phd-scripts/Unpublished/XFEM2/JQuadX.for

293 lines
No EOL
11 KiB
Fortran
Raw Blame History

This file contains ambiguous Unicode characters

This file contains Unicode characters that might be confused with other characters. If you think that this is intentional, you can safely ignore this warning. Use the Escape button to reveal them.

c User subroutine UEL XFEM
subroutine uel(rhs,amatrx,svars,energy,ndofel,nrhs,nsvars,
1 props,nprops,coords,mcrd,nnode,u,du,v,a,jtype,time,dtime,
1 kstep,kinc,jelem,params,ndload,jdltyp,adlmag,predef,npredf,
1 lflags,mlvarx,ddlmag,mdload,pnewdt,jprops,njprop,period)
c
include 'aba_param.inc'
c ABAQUS defined variables:
dimension rhs(mlvarx,*), amatrx(ndofel,ndofel), props(*),
1 svars(nsvars), energy(8), coords(mcrd,nnode), u(ndofel),
1 du(mlvarx,*), v(ndofel), a(ndofel), time(2), params(*),
1 jdltyp(mdload,*), adlmag(mdload,*), ddlmag(mdload,*),
1 predef(2,npredf,nnode), lflags(*), jprops(*)
c
character*256 outdir
integer lenoutdir
integer i,j,k,pss,orderq(3),gint,flag,dimens
integer ncracks,maxncp,nelmx,nnodx,typexe(nnode),ix(nnode)
integer,parameter :: mpg=1650
integer,allocatable:: typex(:,:),ncp(:)
real*8 e, nu
real*8 f(ndofel)
real*8 sg(3,mpg),xypg(2,mpg),xe(8),ye(8),xyc0(2),xycprev(2)
real*8, allocatable:: xyc(:,:,:),dist(:,:),elemgg(:,:)
real*8, allocatable:: batg(:,:),dbatg(:,:),jatg(:)
c
c Read real and integer properties set at the ABAQUS input file
e = props(1)
nu = props(2)
pss = jprops(1)
orderq(1) = jprops(2)
orderq(2) = jprops(3)
orderq(3) = jprops(4)
dimens = jprops(5)
c Read the working directory
call getoutdir(outdir,lenoutdir)
c read number of cracks, max number of crack path points,
c number of enriched elements and enriched nodes.
open(68,file=outdir(1:lenoutdir)//\files\gginfox)
read(68,*) ncracks,maxncp,nelmx,nnodx
close(68)
c Allocate dimensions
allocate (typex(nnodx,2), ncp(ncracks))
allocate (xyc(ncracks,maxncp,2), dist(nnodx,3), elemgg(nelmx,10))
c read coordinates of path points for each crack
open(68,file=outdir(1:lenoutdir)//\files\ggxyc)
do i=1,ncracks
read(68,*) ncp(i)
do j=1,ncp(i)
read(68,*) (xyc(i,j,k),k=1,2)
end do
end do
close(68)
c Read list of enriched nodes, type of enrichment and distances
open(68,file=outdir(1:lenoutdir)//\files\ggnodex)
do i=1,nnodx
read(68,*) (typex(i,j),j=1,2),(dist(i,j),j=2,3)
dist(i,1)=typex(i,1)
end do
close(68)
c read list of enriched elements, type of enrichment and intersection points
open(68,file=outdir(1:lenoutdir)//\files\ggelemx)
do i=1,nelmx
read(68,*) (elemgg(i,j),j=1,10)
end do
close(68)
c call initializing routines for matrix and vectors
call initializem(rhs,ndofel,nrhs)
call initializem(amatrx,ndofel,ndofel)
call initializev(energy,8)
call initializev(svars,nsvars)
c verification of element type (type=12 for enriched element)
if (jtype.eq.12) then
c **************************************
c * 4 node enriched element with *
c * up to 12 dof/node for x-fem *
c **************************************
if (lflags(1).eq.71) then
c coupled thermal-stress, steady state analysis
if (lflags(3).eq.1) then
c Routine that defines the location of integration points according to
c the appropriate subdivision. This enables to know the total number of
c integration points for the current element, stored in gint, and whether
c the element is subdivided for integration (flag=1) or not.
CALL int2d_X(JELEM,NelmX,ElemGG,MCRD,NNODE,COORDS,orderQ,
1 NCracks,maxNCP,NCP,XYC,gint,sg,Xe,Ye,flag,mpg,xypg,
1 XYC0,XYCPrev)
c Allocate dimensions once the total number of integration points gint is known
allocate(batg(3*gint,ndofel),dbatg(3*gint,ndofel),jatg(gint))
call initializem(batg,3*gint,ndofel)
call initializem(dbatg,3*gint,ndofel)
call initializev(jatg,gint)
c Search of the enrichment type for the nodes of the current element.
c The keys to the enrichment types are stored in the element vector TypeXe
call typexelement(outdir,lenoutdir,jelem,nnode,nelmx,ix,typexe)
c element stiffness matrix computation, stored in amatrx
call k_u12(e,nu,amatrx,ndofel,nnode,dimens,mcrd,
coords,pss,nnodx,ix,typexe,dist,xyc0,xycprev,
gint,sg,xe,ye,flag,batg,dbatg,jatg)
c Routine that multiplies AMATRX times U to obtain the force vector F
c at the end of the current increment
call mult_v(amatrx,ndofel,ndofel,u,f,ndofel)
c compute the residual force vector
do i=1,ndofel
rhs(i,1) = rhs(i,1) - f(i)
end do
c Compute stresses at Gauss points for post-processing purposes
c Store them as SVARS for output to the results file (.fil)
call svars_u12(jtype,jelem,svars,nsvars,u,ndofel,batg,
1 dbatg,jatg,gint,mpg,xypg)
end if
end if
end if
return
end
C Element stiffness matrix. Subroutine: K U12
subroutine k_u12(e,nu,amatrx,ndofel,nnode,dimens,mcrd,
1 COORDS,PSS,NnodX,ix,TypeXe,Dist,XYC0,XYCPrev,
1 gint,sg,Xe,Ye,flag,BatG,DBatG,JatG)
implicit none
integer ndofel,nnode,dimens,mcrd,pss,nnodx,gint,flag,pos
integer l,i,j,kk,typexe(nnode),ix(nnode)
real*8 e,nu,dist(nnodx,3),sg(3,*)
real*8 amatrx(ndofel,ndofel),xyc0(2),xycprev(2)
real*8 xe(2*nnode),ye(2*nnode),coords(mcrd,nnode),xl(dimens,nnode)
real*8 xsj(gint),shp(3,4)
real*8 dnf(nnode,2,4),fnode(nnode,4),h,hnode(nnode)
real*8 b(3,ndofel), db(3,ndofel), bt(ndofel,3), d(3,3)
real*8 batg(3*gint,ndofel),dbatg(3*gint,ndofel),jatg(gint)
logical nodetype1,nodetype2
c NOTES:
c Routine shapef2D is called to compute standard shape functions,
c derivatives and jacobian at integration points. This routine outputs:
c shp(3,*) - Shape functions and derivatives at point
c shp(1,i) = dN_i/dx = dN_i/dx1
c shp(2,i) = dN_i/dy = dN_i/dx2
c shp(3,i) = N_i
c xsj - Jacobian determinant at point
c Local coordinates of integration points are passed in sg(1,*), sg(2,*)
c Integration weights are passed in sg(3,*)
c Initialize AMATRX and logical variables
call initializem(amatrx,ndofel,ndofel)
NodeType1=.false.
NodeType2=.false.
c Reduce info passed thru COORDS (3D) to xl (2D)
do i=1,dimens
do j=1,nnode
xl(i,j)=coords(i,j)
end do
end do
c Define constitutive stress-strain elastic matrix
call calc_d(pss,d,e,nu)
c Specify the type of nodal enrichment
do i=1,nnode
if (typexe(i).eq.1) then
nodetype1=.true.
elseif (typexe(i).eq.2) then
nodetype2=.true.
end if
end do
c Numerical integration loop over gint integration points
DO l = 1,gint
c Compute shape functions, derivatives and jacobian at integration point
call shapef2d(sg(1,l),xl,shp,xsj(l),dimens,nnode,ix,.false.)
if (flag.eq.1) then !element is subdivided for integration
xsj(l) = sg(3,l) !the integration weight includes the jacobian
else !element is not subdivided. standard integration
xsj(l) = xsj(l)*sg(3,l)
endif
c Value of the Heaviside function at integration point
c (This call is also used to store the values of H
c at nodes of the element for modified enrichment)
if (nodetype1) then
call heaviside(nnodx,dist,nnode,ix,shp,h,hnode)
endif
c Derivatives of shape functions Ni times enrichment functions Fj at integration point
c (This call is also used to compute the derivatives of shape functions Ni times
c enrichment functions Fj at nodes of the element for modified enrichment)
if (nodetype2) then
call fcracktip(xyc0,xycprev,shp,xe,ye,dnf,fnode)
endif
c STIFFNESS MATRIX COMPUTATION:
c Assembly of element matrix B (denoted as B) at integration point
call initializem(b,3,ndofel)
pos=1
c loop over nodes
do i= 1,nnode
c Contribution to B of derivatives of standard shape functions
B(1,Pos) = shp(1,i)
B(2,Pos+1)= shp(2,i)
B(3,Pos) = shp(2,i)
B(3,Pos+1)= shp(1,i)
c Contribution to B of derivatives of shape functions times Heaviside function
if (typexe(i).eq.1) then
b(1,2+pos) = shp(1,i)*(h-hnode(i))
b(2,3+pos) = shp(2,i)*(h-hnode(i))
b(3,2+pos) = shp(2,i)*(h-hnode(i))
b(3,3+pos) = shp(1,i)*(h-hnode(i))
c Contribution to B of derivatives of shape functions times crack tip functions
elseif(typexe(i).eq.2) then
do kk= 1,4
b(1,2*kk+2+pos)= dnf(i,1,kk)-shp(1,i)*fnode(i,kk)
b(2,2*kk+3+pos)= dnf(i,2,kk)-shp(2,i)*fnode(i,kk)
b(3,2*kk+2+pos)= dnf(i,2,kk)-shp(2,i)*fnode(i,kk)
b(3,2*kk+3+pos)= dnf(i,1,kk)-shp(1,i)*fnode(i,kk)
end do
end if
Pos=Pos+12 !Each node has 12 dof
end do ! i = end loop over element nodes
db=matmul(d,b) ! matrix d*b
bt=transpose(b) ! b transpose
c Integration of BT*D*B
amatrx= amatrx + matmul(bt,db)*xsj(l)
c store information at each integration point for further post-processing
do i=1,3
do j=1,ndofel
batg(3*(l-1)+i,j)=b(i,j)
dbatg(3*(l-1)+i,j)=db(i,j)
end do
end do
jatg(l)=xsj(l)
end do ! l = end loop for each integration point
return
end
c
SUBROUTINE SVARS_U12(JTYPE,JELEM,SVARS,NSVARS,U,Dof,BatG,DBatG,
* JatG,gint,mpg,xypg)
c Calculates and/or stores the following magnitudes at the element integration points,
c storing them in SVARS: strains, stresses, strain energy density, dv/dx, du/dy, jacobian,
c dNi/dx, dNi/dy, global coordinates of integration points.
IMPLICIT NONE
INTEGER i,j,k,NSVARS, Dof, gint, JTYPE,JELEM,mpg
REAL*8 SVARS(NSVARS), U(Dof),BatG(3*gint,Dof),DBatG(3*gint,Dof)
REAL*8 JatG(gint),B(3,Dof),DB(3,Dof),Bdvdx(3,Dof),Bdudy(3,Dof)
REAL*8 EPS(3),SIG(3),W,dvdx(3),dudy(3),JAC,xypg(2,mpg)
c &&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&
39
c &&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&
c First value stored in SVARS is the total number of integration points
c of the enriched element
SVARS(1)=gint
DO i=1,gint
JAC=JatG(i)
DO k=1,3
DO j=1,Dof
B(k,j)=BatG(3*(i-1)+k,j)
Bdvdx(k,j)=B(k,j) ! For computation of dv/dx
Bdudy(k,j)=B(k,j) ! For computation of du/dy
DB(k,j)=DBatG(3*(i-1)+k,j)
END DO
END DO
CALL MULT_V(B,3,Dof,U,EPS,3) ! Compute strains EPS
CALL MULT_V(DB,3,Dof,U,SIG,3) ! Compute stresses SIG
W=0.5d0*(EPS(1)*SIG(1)+EPS(2)*SIG(2)+EPS(3)*SIG(3))
c Computation of dv/dx & du/dy
c Set to zero positions in the 3rd row of B associated with dN/dy
DO j=1,Dof,2
Bdvdx(3,j)=0.0d0
END DO
CALL MULT_V(Bdvdx,3,Dof,U,dvdx,3) !compute dv/dx, stored in dvdx(3)
c Set to zero positions in the 3rd row of B associated with dN/dx
DO j=2,Dof,2
Bdudy(3,j)=0.0d0
END DO
CALL MULT_V(Bdudy,3,Dof,U,dudy,3) !compute du/dy, stored in dudy(3)
c Store in SVARS the following information at integration points
SVARS(1+20*(i-1)+1)=EPS(1)
SVARS(1+20*(i-1)+2)=EPS(2)
SVARS(1+20*(i-1)+3)=EPS(3)
SVARS(1+20*(i-1)+4)=SIG(1)
SVARS(1+20*(i-1)+5)=SIG(2)
SVARS(1+20*(i-1)+6)=SIG(3)
SVARS(1+20*(i-1)+7)=W
SVARS(1+20*(i-1)+8)=dvdx(3)
SVARS(1+20*(i-1)+9)=dudy(3)
SVARS(1+20*(i-1)+10)=JAC ! Jacobian includes integration weight
c Store in SVARS the shape functions derivatives dNi/dx, dNi/dy for external computation
c of dq/dx, dq/dy (used in domain interaction integrals).
c (we take them from the positions associated with the standard dofs)
SVARS(1+20*(i-1)+11)=B(1,1)
SVARS(1+20*(i-1)+12)=B(1,13)
SVARS(1+20*(i-1)+13)=B(1,25)
SVARS(1+20*(i-1)+14)=B(1,37)
SVARS(1+20*(i-1)+15)=B(2,2)
SVARS(1+20*(i-1)+16)=B(2,14)
SVARS(1+20*(i-1)+17)=B(2,26)
SVARS(1+20*(i-1)+18)=B(2,38)
Store in SVARS the global coordinates of integration points
SVARS(1+20*(i-1)+19)=xypg(1,i)
SVARS(1+20*(i-1)+20)=xypg(2,i)
END DO !i loop over all integration points of the element
RETURN
END