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

Go to the source code of this file.

Functions/Subroutines

subroutine networkmpcs (inpc, textpart, ipompc, nodempc, coefmpc, nmpc, nmpc_, mpcfree, nk, ikmpc, ilmpc, labmpc, istep, istat, n, iline, ipol, inl, ipoinp, inp, ipoinpc)
 

Function/Subroutine Documentation

◆ networkmpcs()

subroutine networkmpcs ( character*1, dimension(*)  inpc,
character*132, dimension(16)  textpart,
integer, dimension(*)  ipompc,
integer, dimension(3,*)  nodempc,
real*8, dimension(*)  coefmpc,
integer  nmpc,
integer  nmpc_,
integer  mpcfree,
integer  nk,
integer, dimension(*)  ikmpc,
integer, dimension(*)  ilmpc,
character*20, dimension(*)  labmpc,
integer  istep,
integer  istat,
integer  n,
integer  iline,
integer  ipol,
integer  inl,
integer, dimension(2,*)  ipoinp,
integer, dimension(3,*)  inp,
integer, dimension(0:*)  ipoinpc 
)
22 !
23 ! reading the input deck: *NETWORK MPC
24 !
25  implicit none
26 !
27  character*1 inpc(*)
28  character*13 type
29  character*20 labmpc(*)
30  character*132 textpart(16)
31 !
32  integer ipompc(*),nodempc(3,*),nmpc,nmpc_,mpcfree,istep,istat,
33  & n,i,j,ii,key,nterm,nk,node,ndir,mpcfreeold,ikmpc(*),ilmpc(*),
34  & id,idof,iline,ipol,inl,ipoinp(2,*),inp(3,*),ipoinpc(0:*),m
35 !
36  real*8 coefmpc(*),x
37 !
38  do m=2,n
39  if(textpart(m)(1:5).eq.'TYPE=') then
40  type(1:13)=textpart(m)(6:18)
41  if(textpart(m)(19:19).ne.' ') then
42  write(*,*) '*ERROR reading *NETWORK MPC: type'
43  write(*,*) ' of network mpc is too long'
44  call inputwarning(inpc,ipoinpc,iline,
45  &"*NETWORK MPC%")
46  endif
47  else
48  write(*,*)
49  & '*WARNING reading *NETWORK MPC: parameter not recognized:'
50  write(*,*) ' ',
51  & textpart(m)(1:index(textpart(m),' ')-1)
52  call inputwarning(inpc,ipoinpc,iline,
53  &"*NETWORK MPC%")
54  endif
55  enddo
56 !
57  if(istep.gt.0) then
58  write(*,*)
59  & '*ERROR reading *NETWORK MPC: *NETWORK MPC should be placed'
60  write(*,*) ' before all step definitions'
61  call exit(201)
62  endif
63 !
64  do
65  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
66  & ipoinp,inp,ipoinpc)
67  if((istat.lt.0).or.(key.eq.1)) return
68  read(textpart(1)(1:10),'(i10)',iostat=istat) nterm
69 !
70  nmpc=nmpc+1
71  if(nmpc.gt.nmpc_) then
72  write(*,*) '*ERROR reading *NETWORK MPC: increase nmpc_'
73  call exit(201)
74  endif
75 !
76  labmpc(nmpc)(1:7)='NETWORK'
77  labmpc(nmpc)(8:20)=type(1:13)
78  ipompc(nmpc)=mpcfree
79  ii=0
80 !
81  do
82  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
83  & ipoinp,inp,ipoinpc)
84  if((istat.lt.0).or.(key.eq.1)) then
85  write(*,*)
86  & '*ERROR reading *NETWORK MPC: mpc definition ',nmpc
87  write(*,*) ' is not complete. '
88  call inputerror(inpc,ipoinpc,iline,
89  &"*NETWORK MPC%")
90  call exit(201)
91  endif
92 !
93  do i=1,n/3
94 !
95  read(textpart((i-1)*3+1)(1:10),'(i10)',iostat=istat) node
96  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
97  &"*NETWORK MPC%")
98  if((node.gt.nk).or.(node.le.0)) then
99  write(*,*) '*ERROR reading *NETWORK MPC:'
100  write(*,*) ' node ',node,' is not defined'
101  call exit(201)
102  endif
103 !
104  read(textpart((i-1)*3+2)(1:10),'(i10)',iostat=istat) ndir
105  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
106  &"*NETWORK MPC%")
107  if(ndir.le.6) then
108  elseif(ndir.eq.8) then
109  ndir=4
110  elseif(ndir.eq.11) then
111  ndir=0
112  else
113  write(*,*) '*ERROR reading *NETWORK MPC:'
114  write(*,*) ' direction',ndir,' is not defined'
115  call exit(201)
116  endif
117 !
118  read(textpart((i-1)*3+3)(1:20),'(f20.0)',iostat=istat) x
119  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
120  &"*NETWORK MPC%")
121 !
122  nodempc(1,mpcfree)=node
123  nodempc(2,mpcfree)=ndir
124  coefmpc(mpcfree)=x
125 !
126 ! updating ikmpc and ilmpc
127 !
128  if(ii.eq.0) then
129  idof=8*(node-1)+ndir
130  call nident(ikmpc,idof,nmpc-1,id)
131  if(id.gt.0) then
132  if(ikmpc(id).eq.idof) then
133  write(*,100)
134  & (ikmpc(id))/8+1,ikmpc(id)-8*((ikmpc(id))/8)
135  call exit(201)
136  endif
137  endif
138  do j=nmpc,id+2,-1
139  ikmpc(j)=ikmpc(j-1)
140  ilmpc(j)=ilmpc(j-1)
141  enddo
142  ikmpc(id+1)=idof
143  ilmpc(id+1)=nmpc
144  endif
145 !
146  mpcfreeold=mpcfree
147  mpcfree=nodempc(3,mpcfree)
148  if(mpcfree.eq.0) then
149  write(*,*)
150  & '*ERROR reading *NETWORK MPC: increase memmpc_'
151  call exit(201)
152  endif
153 !
154  ii=ii+1
155  enddo
156 !
157  if(ii.eq.nterm) then
158  nodempc(3,mpcfreeold)=0
159  exit
160  endif
161  enddo
162  enddo
163 !
164  100 format(/,'*ERROR reading *NETWORK MPC: the DOF corresponding to',
165  & /,'node ',i10,' in direction',i1,' is detected on',
166  & /,'the dependent side of two different MPC''s')
167  return
subroutine inputwarning(inpc, ipoinpc, iline, text)
Definition: inputwarning.f:20
subroutine getnewline(inpc, textpart, istat, n, key, iline, ipol, inl, ipoinp, inp, ipoinpc)
Definition: getnewline.f:21
subroutine nident(x, px, n, id)
Definition: nident.f:26
subroutine inputerror(inpc, ipoinpc, iline, text)
Definition: inputerror.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)