User subroutine interface
SUBROUTINE UELMAT(RHS,AMATRX,SVARS,ENERGY,NDOFEL,NRHS,NSVARS,
1 PROPS,NPROPS,COORDS,MCRD,NNODE,U,DU,V,A,JTYPE,TIME,DTIME,
2 KSTEP,KINC,JELEM,PARAMS,NDLOAD,JDLTYP,ADLMAG,PREDEF,NPREDF,
3 LFLAGS,MLVARX,DDLMAG,MDLOAD,PNEWDT,JPROPS,NJPROP,PERIOD,
4 MATERIALLIB)
C
INCLUDE 'ABA_PARAM.INC'
C
DIMENSION RHS(MLVARX,*),AMATRX(NDOFEL,NDOFEL),PROPS(*),
1 SVARS(*),ENERGY(8),COORDS(MCRD,NNODE),U(NDOFEL),
2 DU(MLVARX,*),V(NDOFEL),A(NDOFEL),TIME(2),PARAMS(*),
3 JDLTYP(MDLOAD,*),ADLMAG(MDLOAD,*),DDLMAG(MDLOAD,*),
4 PREDEF(2,NPREDF,NNODE),LFLAGS(*),JPROPS(*)
user coding to define RHS, AMATRX, SVARS, ENERGY, and PNEWDT
RETURN
END
Variables to be defined
These
arrays depend on the value of the LFLAGS
array.
- RHS
-
An array containing the contributions of this element to the right-hand-side
vectors of the overall system of equations. For most nonlinear analysis
procedures, NRHS=1 and
RHS should contain the residual vector. The
exception is the modified Riks static procedure (Static stress analysis),
for which NRHS=2 and the first column in
RHS should contain the residual vector and the
second column should contain the increments of external load on the element.
RHS(K1,K2) is the entry for the
K1th degree of freedom of the element in the
K2th right-hand-side vector.
- AMATRX
-
An array containing the contribution of this element to the Jacobian
(stiffness) or other matrix of the overall system of equations. The particular
matrix required at any time depends on the entries in the
LFLAGS array (see below).
All nonzero entries in AMATRX should be
defined, even if the matrix is symmetric. If you do not specify that the matrix
is unsymmetric when you define the user element,
Abaqus/Standard
will use the symmetric matrix defined by ,
where
is the matrix defined as AMATRX in this
subroutine. If you specify that the matrix is unsymmetric when you define the
user element,
Abaqus/Standard
will use AMATRX directly.
- SVARS
-
An array containing the values of the solution-dependent state variables
associated with this element. The number of such variables is
NSVARS (see below). You define the meaning of
these variables.
For general nonlinear steps this array is passed into
UELMAT containing the values of these variables at the start of
the current increment. They should be updated to be the values at the end of
the increment, unless the procedure during which
UELMAT is being called does not require such an update; this
requirement depends on the entries in the LFLAGS
array (see below). For linear perturbation steps this array is passed into
UELMAT containing the values of these variables in the base
state. They should be returned containing perturbation values if you wish to
output such quantities.
When KINC is equal to zero, the call to
UELMAT is made for zero increment output (see
About Output).
In this case the values returned will be used only for output purposes and are
not updated permanently.
- ENERGY
-
For general nonlinear steps array ENERGY
contains the values of the energy quantities associated with the element. The
values in this array when
UELMAT is called are the element energy quantities at the start
of the current increment. They should be updated to the values at the end of
the current increment. For linear perturbation steps the array is passed into
UELMAT containing the energy in the base state. They should be
returned containing perturbation values if you wish to output such quantities.
The entries in the array are as follows:
ENERGY(1)
|
Kinetic energy.
|
ENERGY(2)
|
Elastic strain energy.
|
ENERGY(3)
|
Creep dissipation.
|
ENERGY(4)
|
Plastic dissipation.
|
ENERGY(5)
|
Viscous dissipation.
|
ENERGY(6)
|
“Artificial strain energy” associated
with such effects as artificial stiffness introduced to control hourglassing or
other singular modes in the element.
|
ENERGY(7)
|
Electrostatic energy.
|
ENERGY(8)
|
Incremental work done by loads applied
within the user element.
|
When KINC is equal to zero, the call to
UELMAT is made for zero increment output (see
About Output).
In this case the energy values returned will be used only for output purposes
and are not updated permanently.
Variables that can be updated
- PNEWDT
-
Ratio of suggested new time increment to the time increment currently being
used (DTIME, see below). This variable allows
you to provide input to the automatic time incrementation algorithms in
Abaqus/Standard
(if automatic time incrementation is chosen). It is useful only during
equilibrium iterations with the normal time incrementation, as indicated by
LFLAGS(3)=1. During a severe discontinuity
iteration (such as contact changes), PNEWDT is
ignored unless CONVERT SDI=YES is specified for this step. The usage of
PNEWDT is discussed below.
PNEWDT is set to a large value before each
call to
UELMAT.
If PNEWDT is redefined to be less than 1.0,
Abaqus/Standard
must abandon the time increment and attempt it again with a smaller time
increment. The suggested new time increment provided to the automatic time
integration algorithms is PNEWDT ×
DTIME, where the
PNEWDT used is the minimum value for all calls
to user subroutines that allow redefinition of
PNEWDT for this iteration.
If PNEWDT is given a value that is greater
than 1.0 for all calls to user subroutines for this iteration and the increment
converges in this iteration,
Abaqus/Standard
may increase the time increment. The suggested new time increment provided to
the automatic time integration algorithms is
PNEWDT × DTIME,
where the PNEWDT used is the minimum value for
all calls to user subroutines for this iteration.
If automatic time incrementation is not selected in the analysis procedure,
values of PNEWDT that are greater than 1.0 will
be ignored and values of PNEWDT that are less
than 1.0 will cause the job to terminate.
Variables passed in for information
- Arrays:
- PROPS
-
A floating point array containing the NPROPS
real property values defined for use with this element.
NPROPS is the user-specified number of real
property values. See
Defining the element properties.
- JPROPS
-
An integer array containing the NJPROP
integer property values defined for use with this element.
NJPROP is the user-specified number of integer
property values. See
Defining the element properties.
- COORDS
-
An array containing the original coordinates of the nodes of the element.
COORDS(K1,K2) is the
K1th coordinate of the
K2th node of the element.
- U, DU, V,
A
-
Arrays containing the current estimates of the basic solution variables
(displacements, rotations, temperatures, depending on the degree of freedom) at
the nodes of the element at the end of the current increment. Values are
provided as follows:
U(K1)
|
Total values of the variables. If this
is a linear perturbation step, it is the value in the base state.
|
DU(K1,KRHS)
|
Incremental values of the variables
for the current increment for right-hand-side
KRHS. If this is an eigenvalue extraction step,
this is the eigenvector magnitude for eigenvector
KRHS. For steady-state dynamics
KRHS
denotes real components of perturbation displacement and
KRHS
denotes imaginary components of perturbation displacement.
|
V(K1)
|
Time rate of change of the variables
(velocities, rates of rotation). Defined for implicit dynamics only
(LFLAGS(1)
11 or 12).
|
A(K1)
|
Accelerations of the variables.
Defined for implicit dynamics only
(LFLAGS(1)
11 or 12).
|
- JDLTYP
-
An array containing the integers used to define distributed load types for
the element. Loads of type Un are identified by the integer value n in
JDLTYP; loads of type UnNU are identified by the negative integer value
in JDLTYP.
JDLTYP(K1,K2) is the identifier of the
K1th distributed load in the
K2th load case. For general nonlinear steps
K2 is always 1.
- ADLMAG
-
For general nonlinear steps ADLMAG(K1,1) is
the total load magnitude of the K1th distributed
load at the end of the current increment for distributed loads of type Un. For distributed loads of type UnNU, the load magnitude is defined in
UELMAT; therefore, the corresponding entries in
ADLMAG are zero. For linear perturbation steps
ADLMAG(K1,1) contains the total load magnitude
of the K1th distributed load of type Un applied in the base state. Base state loading of type UnNU must be dealt with inside
UELMAT. ADLMAG(K1,2),
ADLMAG(K1,3), etc. are currently not used.
- DDLMAG
-
For general nonlinear steps DDLMAG contains
the increments in the magnitudes of the distributed loads that are currently
active on this element for distributed loads of type Un. DDLMAG(K1,1) is the increment of
magnitude of the load for the current time increment. The increment of load
magnitude is needed to compute the external work contribution. For distributed
loads of type UnNU the load magnitude is defined in
UELMAT; therefore, the corresponding entries in
DDLMAG are zero. For linear perturbation steps
DDLMAG(K1,K2) contains the perturbation in the
magnitudes of the distributed loads that are currently active on this element
for distributed loads of type Un. K1 denotes the
K1th perturbation load active on the element.
K2 is always 1, except for steady-state
dynamics, where K2=1 for real loads and
K2=2 for imaginary loads. Perturbation loads of
type UnNU must be dealt with inside
UELMAT.
- PREDEF
-
An array containing the values of predefined field variables, such as
temperature in an uncoupled stress/displacement analysis, at the nodes of the
element (Predefined Fields).
The first index of the array, K1, is either 1
or 2, with 1 indicating the value of the field variable at the end of the
increment and 2 indicating the increment in the field variable. The second
index, K2, indicates the variable: the
temperature corresponds to index 1, and the predefined field variables
correspond to indices 2 and above. In cases where temperature is not defined,
the predefined field variables begin with index 1. The third index,
K3, indicates the local node number on the
element.
PREDEF(K1,1,K3)
|
Temperature.
|
PREDEF(K1,2,K3)
|
First predefined field variable.
|
PREDEF(K1,3,K3)
|
Second predefined field variable.
|
Etc.
|
Any other predefined field variable.
|
PREDEF(K1,K2,K3)
|
Total or incremental value of the
K2th predefined field variable at the
K3th node of the element.
|
PREDEF(1,K2,K3)
|
Values of the variables at the end of
the current increment.
|
PREDEF(2,K2,K3)
|
Incremental values corresponding to
the current time increment.
|
- PARAMS
-
An array containing the parameters associated with the solution procedure.
The entries in this array depend on the solution procedure currently being used
when
UELMAT is called, as indicated by the entries in the
LFLAGS array (see below).
For implicit dynamics (LFLAGS(1) = 11 or 12)
PARAMS contains the integration operator values,
as:
PARAMS(1)
|
|
PARAMS(2)
|
|
PARAMS(3)
|
|
- LFLAGS
-
An array containing the flags that define the current solution procedure and
requirements for element calculations. Detailed requirements for the various
Abaqus/Standard
procedures are defined earlier in this section.
LFLAGS(1)
|
Defines the procedure type. See
Results file output format
for the key used for each procedure.
|
LFLAGS(2)=0
|
Small-displacement analysis.
|
LFLAGS(2)=1
|
Large-displacement analysis (nonlinear
geometric effects included in the step; see
General and perturbation procedures).
|
LFLAGS(3)=1
|
Normal implicit time incrementation
procedure. User subroutine
UELMAT must define the residual vector in
RHS and the Jacobian matrix in
AMATRX.
|
LFLAGS(3)=2
|
Define the current stiffness matrix
(AMATRX)
only.
|
LFLAGS(3)=3
|
Define the current damping matrix
(AMATRX)
only.
|
LFLAGS(3)=4
|
Define the current mass matrix
(AMATRX)
only.
Abaqus/Standard
always requests an initial mass matrix at the start of the analysis.
|
LFLAGS(3)=5
|
Define the current residual or load
vector (RHS)
only.
|
LFLAGS(3)=6
|
Define the current mass matrix and the
residual vector for the initial acceleration calculation (or the calculation of
accelerations after impact).
|
LFLAGS(3)=100
|
Define perturbation quantities for
output.
|
LFLAGS(4)=0
|
The step is a general step.
|
LFLAGS(4)=1
|
The step is a linear perturbation
step.
|
LFLAGS(5)=0
|
The current approximations to
,
etc. were based on Newton corrections.
|
LFLAGS(5)=1
|
The current approximations were found
by extrapolation from the previous increment.
|
- TIME(1)
-
Current value of step time or frequency.
- TIME(2)
-
Current value of total time.
- Scalar parameters:
- DTIME
-
Time increment.
- PERIOD
-
Time period of the current step.
- NDOFEL
-
Number of degrees of freedom in the element.
- MLVARX
-
Dimensioning parameter used when several displacement or right-hand-side
vectors are used.
- NRHS
-
Number of load vectors. NRHS is 1 in most
nonlinear problems: it is 2 for the modified Riks static procedure (Static stress analysis),
and it is greater than 1 in some linear analysis procedures and during
substructure generation.
- NSVARS
-
User-defined number of solution-dependent state variables associated with
the element (Defining the number of solution-dependent variables that must be stored within the element).
- NPROPS
-
User-defined number of real property values associated with the element
(Defining the element properties).
- NJPROP
-
User-defined number of integer property values associated with the element
(Defining the element properties).
- MCRD
-
MCRD is defined as the maximum of the
user-defined maximum number of coordinates needed at any node point (Defining the maximum number of coordinates needed at any nodal point)
and the value of the largest active degree of freedom of the user element that
is less than or equal to 3. For example, if you specify that the maximum number
of coordinates is 1 and the active degrees of freedom of the user element are
2, 3, and 6, MCRD will be 3. If you specify that
the maximum number of coordinates is 2 and the active degrees of freedom of the
user element are 11 and 12, MCRD will be 2.
- NNODE
-
User-defined number of nodes on the element (Defining the number of nodes associated with the element).
- JTYPE
-
Integer defining the element type. This is the user-defined integer value
n in element type Un (Assigning an element type key to a user-defined element).
- KSTEP
-
Current step number.
- KINC
-
Current increment number.
- JELEM
-
User-assigned element number.
- NDLOAD
-
Identification number of the distributed load or flux currently active on
this element.
- MDLOAD
-
Total number of distributed loads and/or fluxes defined on this element.
- NPREDF
-
Number of predefined field variables, including temperature. For user
elements
Abaqus/Standard
uses one value for each field variable per node.
- MATERIALLIB
-
A variable that must be passed to the utility routines performing material
point computations.
UELMAT conventions
The solution variables (displacement, velocity, etc.) are arranged on a
node/degree of freedom basis. The degrees of freedom of the first node are
first, followed by the degrees of freedom of the second node, etc.
Usage with general nonlinear procedures
The values of
(and, in direct-integration dynamic steps,
and )
enter user subroutine
UELMAT as their latest approximations at the end of the time
increment; that is, at time .
The values of
enter the subroutine as their values at the beginning of the time increment;
that is, at time t. It is your responsibility to define
suitable time integration schemes to update .
To ensure accurate, stable integration of internal state variables, you can
control the time incrementation via PNEWDT.
The values of
enter the subroutine as the values of the total load magnitude for the
th
distributed load at the end of the increment. Increments in the load magnitudes
are also available.
In the following descriptions of the user element's requirements, it will be
assumed that LFLAGS(3)=1 unless otherwise
stated.
Static analysis (LFLAGS(1)=1,2)
-
.
-
Automatic convergence checks are applied to the force residuals
corresponding to degrees of freedom 1–7.
-
You must define
AMATRX
and RHS
and update the state variables, .
Direct-integration dynamic analysis (LFLAGS(1)=11,
12)
-
Automatic convergence checks are applied to the force residuals
corresponding to degrees of freedom 1–7.
-
LFLAGS(3)=1: Normal time increment.
Either the Hilber-Hughes-Taylor or the backward Euler time integration scheme
will be used. With
set to zero for the backward Euler, both schemes imply
where
and ;
that is, the highest time derivative of
in
and
is ,
so that
Therefore, you must store
as an internal state vector. If half-increment residual calculations are
required, you must also store
as an internal state vector, where
indicates the time at the beginning of the previous increment. For
,
and
is not needed. You must define
AMATRX
where
and .
RHS
must also be defined and the state variables, ,
updated. Although the value of
given in the dynamic step definition is passed into
UELMAT, the value of
can vary from element to element. For example,
can be set to zero for some elements in the model where numerical dissipation
is not desired.
-
LFLAGS(3)=5: Half-increment residual
()
calculation.
Abaqus/Standard
will adjust the time increment so that
(where
is specified in the dynamic step definition). The half-increment residual is
defined as
where
indicates the time at the beginning of the previous increment
(
is a parameter of the Hilber-Hughes-Taylor time integration operator and will
be set to zero if the backward Euler time integration operator is used). You
must define RHS.
To evaluate
and ,
you must calculate .
These half-increment values will not be saved.
DTIME will still contain
(not ).
The values contained in U,
V, A, and
DU are half-increment values.
-
LFLAGS(3)=4: Velocity jump calculation.
Abaqus/Standard
solves
for ,
so you must define
AMATRX.
-
LFLAGS(3)=6: Initial acceleration
calculation.
Abaqus/Standard
solves
for ,
so you must define
AMATRX
and RHS.
Quasi-static analysis
(LFLAGS(1)=21)
Steady-state heat transfer analysis
(LFLAGS(1)=31)
-
The requirements are identical to those of static analysis, except that
the automatic convergence checks are applied to the heat flux residuals
corresponding to degrees of freedom 11, 12, …
Transient heat transfer analysis
(Δθmax) (LFLAGS(1)=32,
33)
-
Automatic convergence checks are applied to the heat flux residuals
corresponding to degrees of freedom 11, 12, …
-
The backward difference scheme is always used for time integration; that
is,
Abaqus/Standard
assumes that ,
where
and so
always. For degrees of freedom 11, 12, …,
will be compared against the user-prescribed maximum allowable nodal
temperature change in an increment, ,
for controlling the time integration accuracy.
-
You need to define
AMATRX,
where
is the heat capacity matrix and
RHS,
and must update the state variables, .
Usage with linear perturbation procedures
General and perturbation procedures
describes the linear perturbation capabilities in
Abaqus/Standard.
Here, base state values of variables will be denoted by
,
,
etc. Perturbation values will be denoted by ,
,
etc.
Abaqus/Standard
will not call user subroutine
UELMAT for the following procedures: eigenvalue buckling
prediction, response spectrum, transient modal dynamic, steady-state dynamic
(modal and direct), and random response.
Static analysis (LFLAGS(1)=1, 2)
-
Abaqus/Standard
will solve
for ,
where
is the base state stiffness matrix and the perturbation load vector,
,
is a linear function of the perturbation loads, ;
that is, .
-
LFLAGS(3)=1: You must define
AMATRX
and RHS.
-
LFLAGS(3)=100: You must compute
perturbations of the internal variables, ,
and define RHS
for output purposes.
Example: Structural user element with
Abaqus
isotropic linearly elastic material
Both a structural and a heat transfer user element have been created to
demonstrate the usage of subroutine
UELMAT. These user-defined elements are applied in a number of
analyses. The following excerpt illustrates how the linearly elastic isotropic
material available in
Abaqus
can be accessed from user subroutine
UELMAT:
...
USER ELEMENT, TYPE=U1, NODES=4, COORDINATES=2, VAR=16,
INTEGRATION=4, TENSOR=PSTRAIN
1,2
ELEMENT, TYPE=U1, ELSET=SOLID
1, 1,2,3,4
...
UEL PROPERTY, ELSET=SOLID, MATERIAL=MAT
...
MATERIAL, NAME=MAT
ELASTIC
7.00E+010, 0.33
The user element defined above is a 4-node, fully integrated plane strain
element, similar to the
AbaqusCPE4 element.
The next excerpt shows the listing of the user subroutine. Inside the
subroutine, a loop over the integration points is performed. For each
integration point the utility routine MATERIAL_LIB_MECH is called, which returns stress and Jacobian at the integration
point. These quantities are used to compute the right-hand-side vector and the
element Jacobian.
c***********************************************************
subroutine uelmat(rhs,amatrx,svars,energy,ndofel,nrhs,
1 nsvars,props,nprops,coords,mcrd,nnode,u,du,
2 v,a,jtype,time,dtime,kstep,kinc,jelem,params,
3 ndload,jdltyp,adlmag,predef,npredf,lflags,mlvarx,
4 ddlmag,mdload,pnewdt,jprops,njpro,period,
5 materiallib)
c
include 'aba_param.inc'
C
dimension rhs(mlvarx,*), amatrx(ndofel, ndofel), props(*),
1 svars(*), energy(*), coords(mcrd, nnode), u(ndofel),
2 du(mlvarx,*), v(ndofel), a(ndofel), time(2), params(*),
3 jdltyp(mdload,*), adlmag(mdload,*), ddlmag(mdload,*),
4 predef(2, npredf, nnode), lflags(*), jprops(*)
parameter (zero=0.d0, dmone=-1.0d0, one=1.d0, four=4.0d0,
1 fourth=0.25d0,gaussCoord=0.577350269d0)
parameter (ndim=2, ndof=2, nshr=1,nnodemax=4,
1 ntens=4, ninpt=4, nsvint=4)
c
c ndim ... number of spatial dimensions
c ndof ... number of degrees of freedom per node
c nshr ... number of shear stress component
c ntens ... total number of stress tensor components
c (=ndi+nshr)
c ninpt ... number of integration points
c nsvint... number of state variables per integration pt
c (strain)
c
dimension stiff(ndof*nnodemax,ndof*nnodemax),
1 force(ndof*nnodemax), shape(nnodemax), dshape(ndim,nnodemax),
2 xjac(ndim,ndim),xjaci(ndim,ndim), bmat(nnodemax*ndim),
3 statevLocal(nsvint),stress(ntens), ddsdde(ntens, ntens),
4 stran(ntens), dstran(ntens), wght(ninpt)
c
dimension predef_loc(npredf),dpredef_loc(npredf),
1 defGrad(3,3),utmp(3),xdu(3),stiff_p(3,3),force_p(3)
dimension coord24(2,4),coords_ip(3)
data coord24 /dmone, dmone,
2 one, dmone,
3 one, one,
4 dmone, one/
c
data wght /one, one, one, one/
c
c*************************************************************
c
c U1 = first-order, plane strain, full integration
c
c State variables: each integration point has nsvint SDVs
c
c isvinc=(npt-1)*nsvint ... integration point counter
c statev(1+isvinc ) ... strain
c
c*************************************************************
if (lflags(3).eq.4) then
do i=1, ndofel
do j=1, ndofel
amatrx(i,j) = zero
end do
amatrx(i,i) = one
end do
goto 999
end if
c
c PRELIMINARIES
c
pnewdtLocal = pnewdt
if(jtype .ne. 1) then
write(7,*)'Incorrect element type'
call xit
endif
if(nsvars .lt. ninpt*nsvint) then
write(7,*)'Increase the number of SDVs to', ninpt*nsvint
call xit
endif
thickness = 0.1d0
c
c INITIALIZE RHS AND LHS
c
do k1=1, ndof*nnode
rhs(k1, 1)= zero
do k2=1, ndof*nnode
amatrx(k1, k2)= zero
end do
end do
c
c LOOP OVER INTEGRATION POINTS
c
do kintk = 1, ninpt
c
c EVALUATE SHAPE FUNCTIONS AND THEIR DERIVATIVES
c
c determine (g,h)
c
g = coord24(1,kintk)*gaussCoord
h = coord24(2,kintk)*gaussCoord
c
c shape functions
shape(1) = (one - g)*(one - h)/four;
shape(2) = (one + g)*(one - h)/four;
shape(3) = (one + g)*(one + h)/four;
shape(4) = (one - g)*(one + h)/four;
c
c derivative d(Ni)/d(g)
dshape(1,1) = -(one - h)/four;
dshape(1,2) = (one - h)/four;
dshape(1,3) = (one + h)/four;
dshape(1,4) = -(one + h)/four;
c
c derivative d(Ni)/d(h)
dshape(2,1) = -(one - g)/four;
dshape(2,2) = -(one + g)/four;
dshape(2,3) = (one + g)/four;
dshape(2,4) = (one - g)/four;
c
c compute coordinates at the integration point
c
do k1=1, 3
coords_ip(k1) = zero
end do
do k1=1,nnode
do k2=1,mcrd
coords_ip(k2)=coords_ip(k2)+shape(k1)*coords(k2,k1)
end do
end do
c
c INTERPOLATE FIELD VARIABLES
c
if(npredf.gt.0) then
do k1=1,npredf
predef_loc(k1) = zero
dpredef_loc(k1) = zero
do k2=1,nnode
predef_loc(k1) =
& predef_loc(k1)+
& (predef(1,k1,k2)-predef(2,k1,k2))*shape(k2)
dpredef_loc(k1) =
& dpredef_loc(k1)+predef(2,k1,k2)*shape(k2)
end do
end do
end if
c
c FORM B-MATRIX
c
djac = one
c
do i = 1, ndim
do j = 1, ndim
xjac(i,j) = zero
xjaci(i,j) = zero
end do
end do
c
do inod= 1, nnode
do idim = 1, ndim
do jdim = 1, ndim
xjac(jdim,idim) = xjac(jdim,idim) +
1 dshape(jdim,inod)*coords(idim,inod)
end do
end do
end do
djac = xjac(1,1)*xjac(2,2) - xjac(1,2)*xjac(2,1)
if (djac .gt. zero) then
! jacobian is positive - o.k.
xjaci(1,1) = xjac(2,2)/djac
xjaci(2,2) = xjac(1,1)/djac
xjaci(1,2) = -xjac(1,2)/djac
xjaci(2,1) = -xjac(2,1)/djac
else
! negative or zero jacobian
write(7,*)'WARNING: element',jelem,'has neg.
1 Jacobian'
pnewdt = fourth
endif
if (pnewdt .lt. pnewdtLocal) pnewdtLocal = pnewdt
c
do i = 1, nnode*ndim
bmat(i) = zero
end do
do inod = 1, nnode
do ider = 1, ndim
do idim = 1, ndim
irow = idim + (inod - 1)*ndim
bmat(irow) = bmat(irow) +
1 xjaci(idim,ider)*dshape(ider,inod)
end do
end do
end do
c
c CALCULATE INCREMENTAL STRAINS
c
do i = 1, ntens
dstran(i) = zero
end do
!
! set deformation gradient to Identity matrix
do k1=1,3
do k2=1,3
defGrad(k1,k2) = zero
end do
defGrad(k1,k1) = one
end do
c
c COMPUTE INCREMENTAL STRAINS
c
do nodi = 1, nnode
incr_row = (nodi - 1)*ndof
do i = 1, ndof
xdu(i)= du(i + incr_row,1)
utmp(i) = u(i + incr_row)
end do
dNidx = bmat(1 + (nodi-1)*ndim)
dNidy = bmat(2 + (nodi-1)*ndim)
dstran(1) = dstran(1) + dNidx*xdu(1)
dstran(2) = dstran(2) + dNidy*xdu(2)
dstran(4) = dstran(4) +
1 dNidy*xdu(1) +
2 dNidx*xdu(2)
c deformation gradient
defGrad(1,1) = defGrad(1,1) + dNidx*utmp(1)
defGrad(1,2) = defGrad(1,2) + dNidy*utmp(1)
defGrad(2,1) = defGrad(2,1) + dNidx*utmp(2)
defGrad(2,2) = defGrad(2,2) + dNidy*utmp(2)
end do
c
c CALL CONSTITUTIVE ROUTINE
c
isvinc= (kintk-1)*nsvint ! integration point increment
c
c prepare arrays for entry into material routines
c
do i = 1, nsvint
statevLocal(i)=svars(i+isvinc)
end do
c
c state variables
c
!DEC$ NOVECTOR
do k1=1,ntens
stran(k1) = statevLocal(k1)
stress(k1) = zero
end do
c
do i=1, ntens
!DEC$ NOVECTOR
do j=1, ntens
ddsdde(i,j) = zero
end do
ddsdde(i,j) = one
enddo
c
c compute characteristic element length
c
celent = sqrt(djac*dble(ninpt))
dvmat = djac*thickness
c
dvdv0 = one
call material_lib_mech(materiallib,stress,ddsdde,
1 stran,dstran,kintk,dvdv0,dvmat,defGrad,
2 predef_loc,dpredef_loc,npredf,celent,coords_ip)
c
do k1=1,ntens
statevLocal(k1) = stran(k1) + dstran(k1)
end do
c
isvinc= (kintk-1)*nsvint ! integration point increment
c
c update element state variables
c
do i = 1, nsvint
svars(i+isvinc)=statevLocal(i)
end do
c
c form stiffness matrix and internal force vector
c
dNjdx = zero
dNjdy = zero
do i = 1, ndof*nnode
force(i) = zero
do j = 1, ndof*nnode
stiff(j,i) = zero
end do
end do
dvol= wght(kintk)*djac
do nodj = 1, nnode
incr_col = (nodj - 1)*ndof
dNjdx = bmat(1+(nodj-1)*ndim)
dNjdy = bmat(2+(nodj-1)*ndim)
force_p(1) = dNjdx*stress(1) + dNjdy*stress(4)
force_p(2) = dNjdy*stress(2) + dNjdx*stress(4)
do jdof = 1, ndof
jcol = jdof + incr_col
force(jcol) = force(jcol) +
& force_p(jdof)*dvol
end do
do nodi = 1, nnode
incr_row = (nodi -1)*ndof
dNidx = bmat(1+(nodi-1)*ndim)
dNidy = bmat(2+(nodi-1)*ndim)
stiff_p(1,1) = dNidx*ddsdde(1,1)*dNjdx
& + dNidy*ddsdde(4,4)*dNjdy
& + dNidx*ddsdde(1,4)*dNjdy
& + dNidy*ddsdde(4,1)*dNjdx
stiff_p(1,2) = dNidx*ddsdde(1,2)*dNjdy
& + dNidy*ddsdde(4,4)*dNjdx
& + dNidx*ddsdde(1,4)*dNjdx
& + dNidy*ddsdde(4,2)*dNjdy
stiff_p(2,1) = dNidy*ddsdde(2,1)*dNjdx
& + dNidx*ddsdde(4,4)*dNjdy
& + dNidy*ddsdde(2,4)*dNjdy
& + dNidx*ddsdde(4,1)*dNjdx
stiff_p(2,2) = dNidy*ddsdde(2,2)*dNjdy
& + dNidx*ddsdde(4,4)*dNjdx
& + dNidy*ddsdde(2,4)*dNjdx
& + dNidx*ddsdde(4,2)*dNjdy
do jdof = 1, ndof
icol = jdof + incr_col
do idof = 1, ndof
irow = idof + incr_row
stiff(irow,icol) = stiff(irow,icol) +
& stiff_p(idof,jdof)*dvol
end do
end do
end do
end do
c
c assemble rhs and lhs
c
do k1=1, ndof*nnode
rhs(k1, 1) = rhs(k1, 1) - force(k1)
do k2=1, ndof*nnode
amatrx(k1, k2) = amatrx(k1, k2) + stiff(k1,k2)
end do
end do
end do ! end loop on material integration points
pnewdt = pnewdtLocal
c
999 continue
c
return
end
|