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

Go to the source code of this file.

Functions/Subroutines

subroutine submodels (inpc, textpart, set, istartset, iendset, ialset, nset, nset_, nalset, nalset_, nk, istep, istat, n, iline, ipol, inl, ipoinp, inp, ipoinpc, nsubmodel, tieset, tietol, ntie, ntie_, jobnamec, amta, namta, nam, nam_, namtot_)
 

Function/Subroutine Documentation

◆ submodels()

subroutine submodels ( character*1, dimension(*)  inpc,
character*132, dimension(16)  textpart,
character*81, dimension(*)  set,
integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
integer  nset,
integer  nset_,
integer  nalset,
integer  nalset_,
integer  nk,
integer  istep,
integer  istat,
integer  n,
integer  iline,
integer  ipol,
integer  inl,
integer, dimension(2,*)  ipoinp,
integer, dimension(3,*)  inp,
integer, dimension(0:*)  ipoinpc,
integer  nsubmodel,
character*81, dimension(3,*)  tieset,
real*8, dimension(3,*)  tietol,
integer  ntie,
integer  ntie_,
character*132, dimension(*)  jobnamec,
real*8, dimension(2,*)  amta,
integer, dimension(3,*)  namta,
integer  nam,
integer  nam_,
integer  namtot_ 
)
23 !
24 ! reading the input deck: *SUBMODEL
25 !
26  implicit none
27 !
28  character*1 type,inpc(*)
29  character*81 set(*),noset,submset,globset,facialset,
30  & tieset(3,*)
31  character*132 textpart(16),jobnamec(*)
32 !
33  integer nset,nset_,nalset,nalset_,istep,istat,n,key,i,nk,
34  & j,istartset(*),iendset(*),ialset(*),ipos,iline,ipol,inl,
35  & ipoinp(2,*),inp(3,*),l,k,kstart,kend,ipoinpc(0:*),
36  & iglobset,kincrement,nsubmodel,kflag,idummy,nalsetold,ntie,ntie_,
37  & nlength,namta(3,*),nam,nam_,namtot_
38 !
39  real*8 tietol(3,*),amta(2,*)
40 !
41  data kflag /1/
42 !
43  if(istep.gt.0) then
44  write(*,*)
45  & '*ERROR reading *SUBMODEL: *SURFACE should be placed'
46  write(*,*) ' before all step definitions'
47  call exit(201)
48  endif
49 !
50  ntie=ntie+1
51  if(ntie.gt.ntie_) then
52  write(*,*) '*ERROR in ties: increase ntie_'
53  call exit(201)
54  endif
55 !
56  tietol(1,ntie)=-1.d0
57  iglobset=0
58 !
59 ! if no *AMPLITUDE in the input deck, add a fake amplitude
60 ! iamboun and iamload are used to store the global step used
61 ! for interpolation; these fields are not allocated if the
62 ! number of amplitudes is zero
63 !
64  if(nam.eq.0) then
65  nam=1
66  if(nam_.lt.1) then
67  write(*,*) '*ERROR reading *SUBMODEL: increase nam_'
68  call exit(201)
69  endif
70  if(namtot_.lt.2) then
71  write(*,*)
72  & '*ERROR reading *SUBMODEL: increase namtot_'
73  call exit(201)
74  endif
75  amta(1,1)=0.d0
76  amta(2,1)=1.d0
77  amta(1,2)=1.d0
78  amta(2,2)=1.d0
79  namta(1,1)=1
80  namta(2,1)=2
81  namta(3,1)=1
82  endif
83 !
84  type='N'
85 c globset(1:1)=' '
86  nsubmodel=nsubmodel+1
87  submset(1:8)='SUBMODEL'
88  if(nsubmodel.lt.10) then
89  submset(9:10)='00'
90  write(submset(11:11),'(i1)') nsubmodel
91  elseif(nsubmodel.lt.100) then
92  submset(9:9)='0'
93  write(submset(10:11),'(i2)') nsubmodel
94  elseif(nsubmodel.lt.1000) then
95  write(submset(9:11),'(i3)') nsubmodel
96  else
97  write(*,*) '*ERROR reading *SUBMODEL: no more than 999'
98  write(*,*) ' submodels allowed'
99  call exit(201)
100  endif
101  do i=12,81
102  submset(i:i)=' '
103  enddo
104 !
105  do i=2,n
106  if(textpart(i)(1:5).eq.'TYPE=') then
107  if(textpart(i)(6:12).eq.'SURFACE') then
108  type='T'
109  elseif(textpart(i)(6:9).eq.'NODE') then
110  type='N'
111  else
112  write(*,*)
113  & '*ERROR reading *SUBMODEL: unknown type'
114  call exit(201)
115  endif
116  elseif(textpart(i)(1:12).eq.'GLOBALELSET=') then
117  globset(1:80)=textpart(i)(13:92)
118  globset(81:81)=' '
119  ipos=index(globset,' ')
120  globset(ipos:ipos)='E'
121  do iglobset=1,nset
122  if(set(iglobset).eq.globset) exit
123  enddo
124  if(iglobset.gt.nset) then
125  write(*,*)
126  & '*ERROR reading *SUBMODEL: global element set ',
127  & globset(1:ipos-1)
128  write(*,*) ' does not exist'
129  call exit(201)
130  endif
131  do j=istartset(iglobset),iendset(iglobset)
132  if(ialset(j).lt.0) then
133  write(*,*)
134  & '*ERROR reading *SUBMODEL: global element set ',
135  & globset(1:ipos-1)
136  write(*,*) ' was defined using GENERATE;'
137  write(*,*) ' this is not allowed'
138  call exit(201)
139  endif
140  enddo
141  elseif(textpart(i)(1:6).eq.'INPUT=') then
142  jobnamec(4)(1:126)=textpart(i)(7:132)
143  jobnamec(4)(127:132)=' '
144  loop1: do j=1,126
145  if(jobnamec(4)(j:j).eq.'"') then
146  do k=j+1,126
147  if(jobnamec(4)(k:k).eq.'"') then
148  do l=k-1,126
149  jobnamec(4)(l:l)=' '
150  exit loop1
151  enddo
152  endif
153  jobnamec(4)(k-1:k-1)=jobnamec(4)(k:k)
154  enddo
155  jobnamec(4)(126:126)=' '
156  endif
157  enddo loop1
158  else
159  write(*,*)
160  & '*WARNING reading *SUBMODEL: parameter not recognized:'
161  write(*,*) ' ',
162  & textpart(i)(1:index(textpart(i),' ')-1)
163  call inputwarning(inpc,ipoinpc,iline,
164  &"*SUBMODEL%")
165  endif
166  enddo
167 !
168 ! check whether a global elset was defined
169 !
170 c if(globset(1:1).eq.' ') then
171 c write(*,*) '*ERROR reading *SUBMODEL: no global elset'
172 c write(*,*) ' defined'
173 c call exit(201)
174 c endif
175 !
176 c ipos=index(submset,' ')
177 c if(ipos.eq.1) then
178 c write(*,*) '*ERROR reading *SUBMODEL: no name specified'
179 c call exit(201)
180 c endif
181  submset(12:12)=type
182 !
183 ! submodel set (nodal set or element face set)
184 !
185  nset=nset+1
186  if(nset.gt.nset_) then
187  write(*,*) '*ERROR in submsets: increase nset_'
188  call exit(201)
189  endif
190  set(nset)=submset
191  istartset(nset)=nalset+1
192  iendset(nset)=0
193 !
194  if(type.eq.'N') then
195 !
196 ! node surface
197 !
198  do
199  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
200  & ipoinp,inp,ipoinpc)
201  if((istat.lt.0).or.(key.eq.1)) then
202  if(iendset(nset).eq.0) then
203  write(*,*) '*ERROR reading *SUBMODEL: nodal'
204  write(*,*) ' submodel contains no nodes'
205  call exit(201)
206  endif
207  exit
208  endif
209 !
210  do l=1,n
211  if(nalset+1.gt.nalset_) then
212  write(*,*)
213  & '*ERROR reading *SUBMODEL: increase nalset_'
214  call exit(201)
215  endif
216 !
217  read(textpart(l)(1:10),'(i10)',iostat=istat)
218  & ialset(nalset+1)
219  if(istat.gt.0) then
220  noset=textpart(l)(1:80)
221  noset(81:81)=' '
222  ipos=index(noset,' ')
223  noset(ipos:ipos)='N'
224  do i=1,nset
225  if(set(i).eq.noset) then
226  do j=istartset(i),iendset(i)
227  if(ialset(j).gt.0) then
228  nalset=nalset+1
229  if(nalset.gt.nalset_) then
230  write(*,*)
231  & '*ERROR reading *SUBMODEL: increase nalset_'
232  call exit(201)
233  endif
234  ialset(nalset)=ialset(j)
235  else
236  kstart=ialset(nalset-1)
237  kend=ialset(nalset)
238  nalset=nalset-1
239  kincrement=-ialset(j)
240  do k=kstart+kincrement,kend,kincrement
241  nalset=nalset+1
242  if(nalset.gt.nalset_) then
243  write(*,*)
244  & '*ERROR reading *SUBMODEL: increase nalset_'
245  call exit(201)
246  endif
247  ialset(nalset)=k
248  enddo
249  endif
250  enddo
251  iendset(nset)=nalset
252  exit
253  endif
254  enddo
255  if(i.gt.nset-1) then
256  noset(ipos:ipos)=' '
257  write(*,*) '*ERROR reading *SUBMODEL: node set ',
258  & noset
259  write(*,*) ' does not exist'
260  call exit(201)
261  endif
262  else
263  if(ialset(nalset+1).gt.nk) then
264  write(*,*) '*WARNING reading *SUBMODEL: value ',
265  & ialset(nalset+1)
266  write(*,*) ' in set ',set(nset),' > nk'
267  else
268  nalset=nalset+1
269  iendset(nset)=nalset
270  endif
271  endif
272  enddo
273  enddo
274 !
275  else
276 !
277 ! facial submodel interface surface
278 !
279  loop:do
280  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
281  & ipoinp,inp,ipoinpc)
282 !
283 ! check whether any facial surface was defined
284 !
285  if((istat.lt.0).or.(key.eq.1)) then
286  if(iendset(nset).eq.0) then
287  write(*,*) '*ERROR reading *SUBMODEL: facial'
288  write(*,*) ' submodel contains no faces'
289  call exit(201)
290  endif
291  exit
292  endif
293 !
294  if(n.gt.1) then
295  write(*,*) '*ERROR reading *SUBMODEL: no more than'
296  write(*,*) ' one entry per line allowed for'
297  write(*,*) ' a facial submodel interface'
298  call exit(201)
299  endif
300 !
301  facialset(1:80)=textpart(1)(1:80)
302  facialset(81:81)=' '
303  ipos=index(facialset,' ')
304  facialset(ipos:ipos)='T'
305  do i=1,nset
306  if(set(i).eq.facialset) then
307  do j=istartset(i),iendset(i)
308  nalset=nalset+1
309  if(nalset.gt.nalset_) then
310  write(*,*)
311  & '*ERROR in reading *SUBMODEL: increase nalset_'
312  call exit(201)
313  endif
314  ialset(nalset)=ialset(j)
315  enddo
316  iendset(nset)=nalset
317  cycle loop
318  endif
319  enddo
320  enddo loop
321  endif
322 !
323  tieset(1,ntie)(81:81)='S'
324 !
325 ! sorting the local (nodal or element face) set
326 !
327  write(tieset(2,ntie)(1:10),'(i10)') nset
328  tieset(2,ntie)(11:11)=type
329  nlength=iendset(nset)-istartset(nset)+1
330  call isortii(ialset(istartset(nset)),idummy,nlength,kflag)
331 !
332 ! expanding and sorting the global element set
333 !
334  write(tieset(3,ntie)(1:10),'(i10)') iglobset
335 !
336 ! if iglobset equals zero (i.e. no global element set was
337 ! defined), the complete global model is taken as the master model
338 !
339  if(iglobset.eq.0) return
340 !
341 ! expanding the global element set (no "generate" allowed)
342 !
343  nalsetold=nalset+1
344  do j=istartset(iglobset),iendset(iglobset)
345  if(ialset(j).gt.0) then
346  nalset=nalset+1
347  if(nalset.gt.nalset_) then
348  write(*,*)
349  & '*ERROR in surfaces: increase nalset_'
350  call exit(201)
351  endif
352  ialset(nalset)=ialset(j)
353  else
354  kstart=ialset(nalset-1)
355  kend=ialset(nalset)
356  nalset=nalset-1
357  kincrement=-ialset(j)
358  do k=kstart+kincrement,kend,kincrement
359  nalset=nalset+1
360  if(nalset.gt.nalset_) then
361  write(*,*)
362  & '*ERROR in surfaces: increase nalset_'
363  call exit(201)
364  endif
365  ialset(nalset)=k
366  enddo
367  endif
368  enddo
369  istartset(iglobset)=nalsetold
370  iendset(iglobset)=nalset
371 !
372 ! sorting the global element set
373 !
374  nlength=iendset(iglobset)-istartset(iglobset)+1
375  call isortii(ialset(istartset(iglobset)),idummy,nlength,kflag)
376 !
377  return
subroutine inputwarning(inpc, ipoinpc, iline, text)
Definition: inputwarning.f:20
subroutine isortii(ix, iy, n, kflag)
Definition: isortii.f:6
subroutine getnewline(inpc, textpart, istat, n, key, iline, ipol, inl, ipoinp, inp, ipoinpc)
Definition: getnewline.f:21
Hosted by OpenAircraft.com, (Michigan UAV, LLC)