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

Go to the source code of this file.

Functions/Subroutines

subroutine boundaryfs (inpc, textpart, set, istartset, iendset, ialset, nset, nodeboun, ndirboun, xboun, nboun, nboun_, nk, iamboun, amname, nam, ipompc, nodempc, coefmpc, nmpc, nmpc_, mpcfree, trab, ntrans, ikboun, ilboun, ikmpc, ilmpc, nk_, co, labmpc, typeboun, istat, n, iline, ipol, inl, ipoinp, inp, nam_, namtot_, namta, amta, nmethod, iperturb, ipoinpc, vold, mi, xload, sideload, nload, nelemload, lakon, kon, ipkon, ne)
 

Function/Subroutine Documentation

◆ boundaryfs()

subroutine boundaryfs ( 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, dimension(*)  nodeboun,
integer, dimension(*)  ndirboun,
real*8, dimension(*)  xboun,
integer  nboun,
integer  nboun_,
integer  nk,
integer, dimension(*)  iamboun,
character*80, dimension(*)  amname,
integer  nam,
integer, dimension(*)  ipompc,
integer, dimension(3,*)  nodempc,
real*8, dimension(*)  coefmpc,
integer  nmpc,
integer  nmpc_,
integer  mpcfree,
real*8, dimension(7,*)  trab,
integer  ntrans,
integer, dimension(*)  ikboun,
integer, dimension(*)  ilboun,
integer, dimension(*)  ikmpc,
integer, dimension(*)  ilmpc,
integer  nk_,
real*8, dimension(3,*)  co,
character*20, dimension(*)  labmpc,
character*1, dimension(*)  typeboun,
integer  istat,
integer  n,
integer  iline,
integer  ipol,
integer  inl,
integer, dimension(2,*)  ipoinp,
integer, dimension(3,*)  inp,
integer  nam_,
integer  namtot_,
integer, dimension(3,*)  namta,
real*8, dimension(2,*)  amta,
integer  nmethod,
integer  iperturb,
integer, dimension(0:*)  ipoinpc,
real*8, dimension(0:mi(2),*)  vold,
integer, dimension(*)  mi,
real*8, dimension(2,*)  xload,
character*20, dimension(*)  sideload,
integer  nload,
integer, dimension(2,*)  nelemload,
character*8, dimension(*)  lakon,
integer, dimension(*)  kon,
integer, dimension(*)  ipkon,
integer  ne 
)
27 !
28 ! reading the input deck: *BOUNDARYF
29 !
30 ! (boundary conditions for CFD-calculations)
31 !
32  implicit none
33 !
34  logical user,fixed,surface
35 !
36  character*1 typeboun(*),type,inpc(*)
37  character*8 lakon(*)
38  character*20 labmpc(*),sideload(*)
39  character*80 amname(*),amplitude
40  character*81 set(*),elset
41  character*132 textpart(16)
42 !
43  integer istartset(*),iendset(*),ialset(*),nodeboun(*),
44  & ndirboun(*),iface,nload,nelemload(2,*),kon(*),ipkon(*),
45  & nset,nboun,nboun_,istat,n,i,j,k,l,ibounstart,ibounend,
46  & key,nk,iamboun(*),nam,iamplitude,ipompc(*),nodempc(3,*),
47  & nmpc,nmpc_,mpcfree,ikboun(*),ilboun(*),ikmpc(*),
48  & ilmpc(*),ntrans,nk_,ipos,m,ne,
49  & iline,ipol,inl,ipoinp(2,*),inp(3,*),nam_,namtot,namtot_,
50  & namta(3,*),idelay,nmethod,iperturb,ipoinpc(0:*),
51  & mi(*)
52 !
53  real*8 xboun(*),bounval,coefmpc(*),trab(7,*),co(3,*),amta(2,*),
54  & vold(0:mi(2),*),xload(2,*)
55 !
56  type='F'
57  iamplitude=0
58  idelay=0
59  user=.false.
60  fixed=.false.
61  surface=.false.
62 !
63  do i=2,n
64  if(textpart(i)(1:10).eq.'AMPLITUDE=') then
65  read(textpart(i)(11:90),'(a80)') amplitude
66  do j=nam,1,-1
67  if(amname(j).eq.amplitude) then
68  iamplitude=j
69  exit
70  endif
71  enddo
72  if(j.eq.0) then
73  write(*,*)
74  & '*ERROR reading *BOUNDARYF: nonexistent amplitude'
75  write(*,*) ' '
76  call inputerror(inpc,ipoinpc,iline,
77  &"*BOUNDARYF%")
78  call exit(201)
79  endif
80  iamplitude=j
81  elseif(textpart(i)(1:10).eq.'TIMEDELAY=') THEN
82  if(idelay.ne.0) then
83  write(*,*)'*ERROR reading *BOUNDARYF: the parameter TIME'
84  write(*,*) ' DELAY is used twice in the same'
85  write(*,*) ' keyword; '
86  call inputerror(inpc,ipoinpc,iline,
87  &"*BOUNDARYF%")
88  call exit(201)
89  else
90  idelay=1
91  endif
92  nam=nam+1
93  if(nam.gt.nam_) then
94  write(*,*) '*ERROR reading *BOUNDARYF: increase nam_'
95  call exit(201)
96  endif
97  amname(nam)='
98  & '
99  if(iamplitude.eq.0) then
100  write(*,*)'*ERROR reading *BOUNDARYF: time delay must be'
101  write(*,*) ' preceded by the amplitude parameter'
102  call exit(201)
103  endif
104  namta(3,nam)=sign(iamplitude,namta(3,iamplitude))
105  iamplitude=nam
106  if(nam.eq.1) then
107  namtot=0
108  else
109  namtot=namta(2,nam-1)
110  endif
111  namtot=namtot+1
112  if(namtot.gt.namtot_) then
113  write(*,*) '*ERROR boundaries: increase namtot_'
114  call exit(201)
115  endif
116  namta(1,nam)=namtot
117  namta(2,nam)=namtot
118  read(textpart(i)(11:30),'(f20.0)',iostat=istat)
119  & amta(1,namtot)
120  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
121  &"*BOUNDARYF%")
122  elseif(textpart(i)(1:4).eq.'USER') then
123  user=.true.
124  else
125  write(*,*)
126  & '*WARNING reading *BOUNDARYF: parameter not recognized:'
127  write(*,*) ' ',
128  & textpart(i)(1:index(textpart(i),' ')-1)
129  call inputwarning(inpc,ipoinpc,iline,
130  &"*BOUNDARYF%")
131  endif
132  enddo
133 !
134  if(user.and.(iamplitude.ne.0)) then
135  write(*,*) '*WARNING: no amplitude definition is allowed'
136  write(*,*) ' for temperatures defined by a'
137  write(*,*) ' user routine'
138  iamplitude=0
139  endif
140 !
141  do
142  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
143  & ipoinp,inp,ipoinpc)
144  if((istat.lt.0).or.(key.eq.1)) return
145 !
146  read(textpart(3)(1:10),'(i10)',iostat=istat) ibounstart
147  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
148  &"*BOUNDARYF%")
149 !
150  if(textpart(4)(1:1).eq.' ') then
151  ibounend=ibounstart
152  else
153  read(textpart(4)(1:10),'(i10)',iostat=istat) ibounend
154  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
155  &"*BOUNDARYF%")
156  endif
157 !
158  if(textpart(5)(1:1).eq.' ') then
159  bounval=0.d0
160  else
161  read(textpart(5)(1:20),'(f20.0)',iostat=istat) bounval
162  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
163  &"*BOUNDARYF%")
164  endif
165 !
166 ! dummy boundary condition consisting of the first primes
167 !
168  if(user) bounval=1.2357111317d0
169 !
170  read(textpart(1)(1:10),'(i10)',iostat=istat) l
171  if(istat.eq.0) then
172  if((l.gt.ne).or.(l.le.0)) then
173  write(*,*) '*ERROR reading *BOUNDARYF:'
174  write(*,*) ' element ',l,' is not defined'
175  call exit(201)
176  endif
177  read(textpart(2)(2:2),'(i1)',iostat=istat) iface
178  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
179  &"*BOUNDARYF%")
180  l=10*l+iface
181  call bounaddf(l,ibounstart,ibounend,bounval,
182  & nodeboun,ndirboun,xboun,nboun,nboun_,
183  & iamboun,iamplitude,nam,ipompc,nodempc,
184  & coefmpc,nmpc,nmpc_,mpcfree,trab,
185  & ntrans,ikboun,ilboun,ikmpc,ilmpc,co,nk,nk_,labmpc,
186  & type,typeboun,nmethod,iperturb,vold,mi,
187  & nelemload,sideload,xload,nload,lakon,ipkon,kon)
188  else
189  read(textpart(1)(1:80),'(a80)',iostat=istat) elset
190  elset(81:81)=' '
191  ipos=index(elset,' ')
192  elset(ipos:ipos)='E'
193  do i=1,nset
194  if(set(i).eq.elset) exit
195  enddo
196  if(i.gt.nset) then
197 !
198 ! check for facial surface
199 !
200  surface=.true.
201  elset(ipos:ipos)='T'
202  do i=1,nset
203  if(set(i).eq.elset) exit
204  enddo
205  if(i.gt.nset) then
206  elset(ipos:ipos)=' '
207  write(*,*) '*ERROR reading *BOUNDARYF: surface ',elset
208  write(*,*) ' has not yet been defined. '
209  call inputerror(inpc,ipoinpc,iline,
210  &"*BOUNDARYF%")
211  call exit(201)
212  endif
213  endif
214  read(textpart(2)(2:2),'(i1)',iostat=istat) iface
215  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
216  &"*BOUNDARYF%")
217  do j=istartset(i),iendset(i)
218  if(ialset(j).gt.0) then
219  k=ialset(j)
220  if(.not.surface) k=10*k+iface
221  call bounaddf(k,ibounstart,ibounend,bounval,
222  & nodeboun,ndirboun,xboun,nboun,nboun_,
223  & iamboun,iamplitude,nam,ipompc,nodempc,
224  & coefmpc,nmpc,nmpc_,mpcfree,trab,
225  & ntrans,ikboun,ilboun,ikmpc,ilmpc,co,nk,nk_,labmpc,
226  & type,typeboun,nmethod,iperturb,vold,mi,
227  & nelemload,sideload,xload,nload,lakon,ipkon,kon)
228  else
229  m=ialset(j-2)
230  do
231  m=m-ialset(j)
232  if(m.ge.ialset(j-1)) exit
233  k=10*m+iface
234  call bounaddf(k,ibounstart,ibounend,bounval,
235  & nodeboun,ndirboun,xboun,nboun,nboun_,
236  & iamboun,iamplitude,nam,ipompc,nodempc,
237  & coefmpc,nmpc,nmpc_,mpcfree,trab,
238  & ntrans,ikboun,ilboun,ikmpc,ilmpc,co,nk,nk_,
239  & labmpc,type,typeboun,nmethod,iperturb,
240  & vold,mi,
241  & nelemload,sideload,xload,nload,lakon,ipkon,kon)
242  enddo
243  endif
244  enddo
245  endif
246  enddo
247 !
248  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 inputerror(inpc, ipoinpc, iline, text)
Definition: inputerror.f:20
subroutine bounaddf(iface, is, ie, val, nodeboun, ndirboun, xboun, nboun, nboun_, iamboun, iamplitude, nam, ipompc, nodempc, coefmpc, nmpc, nmpc_, mpcfree, trab, ntrans, ikboun, ilboun, ikmpc, ilmpc, co, nk, nk_, labmpc, type, typeboun, nmethod, iperturb, vold, mi, nelemload, sideload, xload, nload, lakon, ipkon, kon)
Definition: bounaddf.f:25
Hosted by OpenAircraft.com, (Michigan UAV, LLC)