1!
2!     CalculiX - A 3-dimensional finite element program
3!     Copyright (C) 1998-2021 Guido Dhondt
4!
5!     This program is free software; you can redistribute it and/or
6!     modify it under the terms of the GNU General Public License as
7!     published by the Free Software Foundation(version 2);
8!
9!
10!     This program is distributed in the hope that it will be useful,
11!     but WITHOUT ANY WARRANTY; without even the implied warranty of
12!     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13!     GNU General Public License for more details.
14!
15!     You should have received a copy of the GNU General Public License
16!     along with this program; if not, write to the Free Software
17!     Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
18!
19      subroutine e_c3d_u(co,kon,lakonl,p1,p2,omx,bodyfx,nbody,s,sm,
20     &     ff,nelem,nmethod,elcon,nelcon,rhcon,nrhcon,alcon,nalcon,
21     &     alzero,ielmat,ielorien,norien,orab,ntmat_,
22     &     t0,t1,ithermal,vold,iperturb,nelemload,
23     &     sideload,xload,nload,idist,sti,stx,iexpl,plicon,
24     &     nplicon,plkcon,nplkcon,xstiff,npmat_,dtime,
25     &     matname,mi,ncmat_,mass,stiffness,buckling,rhsi,intscheme,
26     &     ttime,time,istep,iinc,coriolis,xloadold,reltime,
27     &     ipompc,nodempc,coefmpc,nmpc,ikmpc,ilmpc,veold,
28     &     ne0,ipkon,thicke,
29     &     integerglob,doubleglob,tieset,istartset,iendset,ialset,ntie,
30     &     nasym,ielprop,prop)
31!
32!     computation of the element matrix and rhs for the user element
33!     of type u1
34!
35!     U1: This is a beam type element. Reference:
36!     Yunhua Luo, An Efficient 3D Timoshenko Beam Element with
37!     Consistent Shape Functions, Adv. Theor. Appl. Mech., Vol. 1,
38!     2008, no. 3, 95-106
39!
40!     special case for which the beam axis goes through the
41!     center of gravity of the cross section and the 1-direction
42!     corresponds with a principal axis
43!
44!     US45: Flat shell 4 nodes 5 dof + (artificial rotz stiff)
45!
46      implicit none
47!
48      integer mass(*),stiffness,buckling,rhsi,coriolis
49!
50      character*8 lakonl
51      character*20 sideload(*)
52      character*80 matname(*),amat
53      character*81 tieset(3,*)
54!
55      integer konl(26),nelemload(2,*),nbody,nelem,mi(*),kon(*),
56     &     ielprop(*),null,index,mattyp,ithermal(*),iperturb(*),nload,
57     &     idist,
58     &     i,j,i1,nmethod,kk,nelcon(2,*),nrhcon(*),nalcon(2,*),
59     &     ielmat(mi(3),*),ielorien(mi(3),*),ipkon(*),indexe,
60     &     ntmat_,nope,norien,ihyper,iexpl,kode,imat,iorien,istiff,
61     &     ncmat_,intscheme,istep,iinc,iflag,ipompc(*),nodempc(3,*),
62     &     nmpc,ikmpc(*),ilmpc(*),ne0,ndof,istartset(*),iendset(*),
63     &     ialset(*),ntie,integerglob(*),nasym,nplicon(0:ntmat_,*),
64     &     nplkcon(0:ntmat_,*),npmat_
65!
66      real*8 co(3,*),xl(3,26),veold(0:mi(2),*),rho,s(60,60),bodyfx(3),
67     &     ff(60),elconloc(21),coords(3),p1(3),elcon(0:ncmat_,ntmat_,*),
68     &     p2(3),eth(6),rhcon(0:1,ntmat_,*),reltime,prop(*),
69     &     alcon(0:6,ntmat_,*),alzero(*),orab(7,*),t0(*),t1(*),
70     &     xloadold(2,*),vold(0:mi(2),*),xload(2,*),omx,e,un,um,tt,
71     &     sm(60,60),sti(6,mi(1),*),stx(6,mi(1),*),t0l,t1l,coefmpc(*),
72     &     elas(21),thicke(mi(3),*),doubleglob(*),dl,
73     &     plicon(0:2*npmat_,ntmat_,*),plkcon(0:2*npmat_,ntmat_,*),
74     &     xstiff(27,mi(1),*),plconloc(802),dtime,ttime,time,
75     &     a,xi11,xi12,xi22,xk,e1(3),offset1,offset2,y1,y2,y3,z1,z2,z3
76!
77      if(lakonl(2:2).eq.'1') then
78        call e_c3d_u1(co,kon,lakonl,p1,p2,omx,bodyfx,nbody,s,sm,
79     &       ff,nelem,nmethod,elcon,nelcon,rhcon,nrhcon,alcon,nalcon,
80     &       alzero,ielmat,ielorien,norien,orab,ntmat_,
81     &       t0,t1,ithermal,vold,iperturb,nelemload,
82     &       sideload,xload,nload,idist,sti,stx,iexpl,plicon,
83     &       nplicon,plkcon,nplkcon,xstiff,npmat_,dtime,
84     &       matname,mi,ncmat_,mass,stiffness,buckling,rhsi,intscheme,
85     &       ttime,time,istep,iinc,coriolis,xloadold,reltime,
86     &       ipompc,nodempc,coefmpc,nmpc,ikmpc,ilmpc,veold,
87     &       ne0,ipkon,thicke,integerglob,doubleglob,tieset,istartset,
88     &       iendset,ialset,ntie,nasym,ielprop,prop)
89!
90      elseif(lakonl(2:4).eq.'S45') then
91        call e_c3d_us45(co,kon,lakonl,p1,p2,omx,bodyfx,nbody,s,sm,
92     &       ff,nelem,nmethod,elcon,nelcon,rhcon,nrhcon,alcon,nalcon,
93     &       alzero,ielmat,ielorien,norien,orab,ntmat_,
94     &       t0,t1,ithermal,vold,iperturb,nelemload,
95     &       sideload,xload,nload,idist,sti,stx,iexpl,plicon,
96     &       nplicon,plkcon,nplkcon,xstiff,npmat_,dtime,
97     &       matname,mi,ncmat_,mass,stiffness,buckling,rhsi,intscheme,
98     &       ttime,time,istep,iinc,coriolis,xloadold,reltime,
99     &       ipompc,nodempc,coefmpc,nmpc,ikmpc,ilmpc,veold,
100     &       ne0,ipkon,thicke,integerglob,doubleglob,tieset,istartset,
101     &       iendset,ialset,ntie,nasym,ielprop,prop)
102!
103      elseif(lakonl(2:3).eq.'S3') then
104        call e_c3d_us3(co,kon,lakonl,p1,p2,omx,bodyfx,nbody,s,sm,
105     &       ff,nelem,nmethod,elcon,nelcon,rhcon,nrhcon,alcon,nalcon,
106     &       alzero,ielmat,ielorien,norien,orab,ntmat_,
107     &       t0,t1,ithermal,vold,iperturb,nelemload,
108     &       sideload,xload,nload,idist,sti,stx,iexpl,plicon,
109     &       nplicon,plkcon,nplkcon,xstiff,npmat_,dtime,
110     &       matname,mi,ncmat_,mass,stiffness,buckling,rhsi,intscheme,
111     &       ttime,time,istep,iinc,coriolis,xloadold,reltime,
112     &       ipompc,nodempc,coefmpc,nmpc,ikmpc,ilmpc,veold,
113     &       ne0,ipkon,thicke,integerglob,doubleglob,tieset,istartset,
114     &       iendset,ialset,ntie,nasym,ielprop,prop)
115      else
116        write(*,*) '*ERROR in e_c3d_u.f: user element'
117        write(*,*) '       ',lakonl(1:5),' is not defined'
118        call exit(201)
119      endif
120!
121      return
122      end
123