added interface routine for Marc2008r1
this includes two files: concom_f90 and creeps_f90 so you either create links of that name to the new files concom2008r1 and creeps2008r1 or rename copies of those files accordingly
This commit is contained in:
parent
fa7cf61fd9
commit
1a5a68c87a
|
@ -0,0 +1,186 @@
|
||||||
|
! reformated to free format
|
||||||
|
!***********************************************************************
|
||||||
|
!
|
||||||
|
! File: concom.cmn
|
||||||
|
!
|
||||||
|
! MSC.Marc include file
|
||||||
|
!
|
||||||
|
integer(pInt) &
|
||||||
|
iacous, iasmbl, iautth, ibear, icompl, iconj, icreep, ideva, idyn, idynt,&
|
||||||
|
ielas, ielcma, ielect, iform, ifour, iharm, ihcps, iheat, iheatt, ihresp,&
|
||||||
|
ijoule, ilem, ilnmom, iloren, inc, incext, incsub, ipass, iplres, ipois,&
|
||||||
|
ipoist, irpflo, ismall, ismalt, isoil, ispect, ispnow, istore, iswep, ithcrp,&
|
||||||
|
itherm, iupblg, iupdat, jacflg, jel, jparks, largst, lfond, loadup, loaduq,&
|
||||||
|
lodcor, lovl, lsub, magnet, ncycle, newtnt, newton, noshr, linear, ivscpl,&
|
||||||
|
icrpim, iradrt, ipshft, itshr, iangin, iupmdr, iconjf, jincfl, jpermg, jhour,&
|
||||||
|
isolvr, jritz, jtable, jshell, jdoubl, jform, jcentr, imini, kautth, iautof,&
|
||||||
|
ibukty, iassum, icnstd, icnstt, kmakmas, imethvp,iradrte,iradrtp, iupdate,iupdatp,&
|
||||||
|
ncycnt, marmen ,idynme, ihavca, ispf, kmini, imixed, largtt, kdoela, iautofg,&
|
||||||
|
ipshftp,idntrc, ipore, jtablm, jtablc, isnecma,itrnspo,imsdif, jtrnspo,mcnear,&
|
||||||
|
imech, imecht, ielcmat, ielectt,magnett, imsdift,noplas, jtabls, jactch, jtablth,&
|
||||||
|
kgmsto ,jpzo, ifricsh, iremkin,iremfor, ishearp,jspf, machining, jlshell,icompsol,&
|
||||||
|
iupblgfo,jcondir,nstcrp, nactive,ipassref, nstspnt,ibeart,icheckmpc, noline, icuring,&
|
||||||
|
ishrink,ioffsflg,isetoff, ioffsetm,iharmt, inc_incdat,iautspc,ibrake, icbush ,istream_input,&
|
||||||
|
iprsinp,ivlsinp,ifirst_time,ipin_m,jgnstr_glb,imarc_return,iqvcinp,nqvceid,istpnx,imicro1,&
|
||||||
|
iaxisymm,jbreakglue,iglstif,jfastasm,iwear, iwearcf, imixmeth,ielcmadyn,idinout,igena_meth
|
||||||
|
integer(pInt) num_concom
|
||||||
|
parameter(num_concom=219)
|
||||||
|
common/marc_concom/&
|
||||||
|
iacous, iasmbl, iautth, ibear, icompl, iconj, icreep, ideva(50), idyn, idynt,&
|
||||||
|
ielas, ielcma, ielect, iform, ifour, iharm, ihcps, iheat, iheatt, ihresp,&
|
||||||
|
ijoule, ilem, ilnmom, iloren, inc, incext, incsub, ipass, iplres, ipois,&
|
||||||
|
ipoist, irpflo, ismall, ismalt, isoil, ispect, ispnow, istore, iswep, ithcrp,&
|
||||||
|
itherm, iupblg, iupdat, jacflg, jel, jparks, largst, lfond, loadup, loaduq,&
|
||||||
|
lodcor, lovl, lsub, magnet, ncycle, newtnt, newton, noshr, linear, ivscpl,&
|
||||||
|
icrpim, iradrt, ipshft, itshr, iangin, iupmdr, iconjf, jincfl, jpermg, jhour,&
|
||||||
|
isolvr, jritz, jtable, jshell, jdoubl, jform, jcentr, imini, kautth, iautof,&
|
||||||
|
ibukty, iassum, icnstd, icnstt, kmakmas, imethvp,iradrte,iradrtp, iupdate,iupdatp,&
|
||||||
|
ncycnt, marmen ,idynme, ihavca, ispf, kmini, imixed, largtt, kdoela, iautofg,&
|
||||||
|
ipshftp,idntrc, ipore, jtablm, jtablc, isnecma,itrnspo,imsdif, jtrnspo,mcnear,&
|
||||||
|
imech, imecht, ielcmat, ielectt,magnett, imsdift,noplas, jtabls, jactch, jtablth,&
|
||||||
|
kgmsto ,jpzo, ifricsh, iremkin,iremfor, ishearp,jspf, machining, jlshell,icompsol,&
|
||||||
|
iupblgfo,jcondir,nstcrp, nactive,ipassref, nstspnt,ibeart,icheckmpc, noline, icuring,&
|
||||||
|
ishrink,ioffsflg,isetoff, ioffsetm,iharmt, inc_incdat,iautspc,ibrake, icbush ,istream_input,&
|
||||||
|
iprsinp,ivlsinp,ifirst_time,ipin_m,jgnstr_glb,imarc_return,iqvcinp,nqvceid,istpnx,imicro1,&
|
||||||
|
iaxisymm,jbreakglue,iglstif,jfastasm,iwear, iwearcf, imixmeth, ielcmadyn,idinout,igena_meth
|
||||||
|
!
|
||||||
|
! comments of variables:
|
||||||
|
!
|
||||||
|
! ideva(50) - debug print out flag
|
||||||
|
! 1 print element stiffness matrices, mass matrix
|
||||||
|
! 2 output matrices used in tying
|
||||||
|
! 3 force the solution of a nonpositive definite matrix
|
||||||
|
! 4 print info of connections to each node
|
||||||
|
! 5 info of gap convergence, internal heat generated, contact
|
||||||
|
! touching and separation
|
||||||
|
! 6 nodal value array during rezoning
|
||||||
|
! 7 tying info in CONRAD GAP option, fluid element numbers in
|
||||||
|
! CHANNEL option
|
||||||
|
! 8 output incremental displacements in local coord. system
|
||||||
|
! 9 latent heat output
|
||||||
|
! 10 stress-strain in local coord. system
|
||||||
|
! 11 additional info on interlaminar stress
|
||||||
|
! 12 output right hand side and solution vector
|
||||||
|
! 13 info of CPU resources used and memory available on NT
|
||||||
|
! 14 info of mesh adaption process, 2D outline information
|
||||||
|
! info of penetration checking for remeshing
|
||||||
|
! save .fem files after afmesh3d meshing
|
||||||
|
! 15 surface energy balance flag
|
||||||
|
! 16 print info regarding pyrolysis
|
||||||
|
! 17 print info of "streamline topology"
|
||||||
|
! 18 print mesh data changes after remeshing
|
||||||
|
! 19 print material flow stress data read in from *.mat file
|
||||||
|
! if unit flag is on, print out flow stress after conversion
|
||||||
|
! 20 print information on table input
|
||||||
|
! 21 print out information regarding kinematic boundary conditions
|
||||||
|
! 22 print out information regarding dist loads, point loads, film
|
||||||
|
! and foundations
|
||||||
|
! 23 print out information about automatic domain decomposition
|
||||||
|
! 24 print out iteration information in SuperForm status report file
|
||||||
|
! 25 print out information for ablation
|
||||||
|
! 26 print out information for films - Table input
|
||||||
|
! 27 print out the tying forces
|
||||||
|
! 28 print out for CASI solver, convection,
|
||||||
|
! 29 DDM single file debug printout
|
||||||
|
! 30 print out cavity debug info
|
||||||
|
! 31 print out welding related info
|
||||||
|
! 32 prints categorized DDM memory usage
|
||||||
|
! 33 print out the cutting info regarding machining feature
|
||||||
|
! 34 print out the list of quantities which can be defined via a table
|
||||||
|
! and for each quantity the supported independent variables
|
||||||
|
! 35 print out detailed coupling region info
|
||||||
|
! 36 print out solver debug info level 1 (Least Detailed)
|
||||||
|
! 37 print out solver debug info level 1 (Medium Detailed)
|
||||||
|
! 38 print out solver debug info level 1 (Very Detailed)
|
||||||
|
! 39 print detailed memory allocation info
|
||||||
|
! 40 print out marc-adams debug info
|
||||||
|
! 41 output rezone mapping post file for debugging
|
||||||
|
! 42 output post file after calling oprofos() for debugging
|
||||||
|
! 43 debug printout for vcct
|
||||||
|
! 44 debug printout for progressive failure
|
||||||
|
! 45 print out automatically generated midside node coordinates (arecrd)
|
||||||
|
! 46 print out message about routine and location, where the ibort is raised (ibort_inc)
|
||||||
|
! 47 print out summary message of element variables on a
|
||||||
|
! group-basis after all the automatic changes have been
|
||||||
|
! made (em_ellibp)
|
||||||
|
! 48 Automatically generate check results based on max and min vals.
|
||||||
|
! These vals are stored in the checkr file, which is inserted
|
||||||
|
! into the *dat file by the generate_check_results script from /marc/tools
|
||||||
|
! 49 Automatically generate check results based on the real calculated values
|
||||||
|
! at the sppecified check result locations.
|
||||||
|
! These vals are stored in the checkr file, which is inserted
|
||||||
|
! into the *dat file by the update_check_results script from /marc/tools
|
||||||
|
!
|
||||||
|
!
|
||||||
|
! jactch = 1 or 2 if elements are activated or deactivated
|
||||||
|
! = 3 if elements are adaptively remeshed or rezoned
|
||||||
|
! = 0 normally / reset to 0 when assembly is done
|
||||||
|
! ifricsh = 0 call to fricsh in otest not needed
|
||||||
|
! = 1 call to fricsh (nodal friction) in otest needed
|
||||||
|
! iremkin = 0 remove deactivated kinematic boundary conditions
|
||||||
|
! immediately - only in new input format (this is default)
|
||||||
|
! = 1 remove deactivated kinematic boundary conditions
|
||||||
|
! gradually - only in new input format
|
||||||
|
! iremfor = 0 remove force boundary conditions immediately -
|
||||||
|
! only in new input format (this is default)
|
||||||
|
! = 1 remove force boundary conditions gradually -
|
||||||
|
! only in new input format (this is default)
|
||||||
|
! ishearp set to 1 if shear panel elements are present in the model
|
||||||
|
!
|
||||||
|
! jspf = 0 not in spf loadcase
|
||||||
|
! > 0 in spf loadcase (jspf=1 during first increment)
|
||||||
|
! machining = 1 if the metal cutting feature is used, for memory allocation purpose
|
||||||
|
! = 0 (default) if no metal cutting feature required
|
||||||
|
!
|
||||||
|
! jlshell = 1 if there is a shell element in the mesh
|
||||||
|
! icompsol = 1 if there is a composite solid element in the mesh
|
||||||
|
! iupblgfo = 1 if follower force for point loads
|
||||||
|
! jcondir = 1 if contact priority option is used
|
||||||
|
! nstcrp = 0 (default) steady state creep flag (undocumented feature.
|
||||||
|
! if not 0, turns off special ncycle = 0 code in radial.f)
|
||||||
|
! nactive = number of active passes, if =1 then it's not a coupled analysis
|
||||||
|
! ipassref = reference ipass, if not in a multiphysics pass ipass=ipassref
|
||||||
|
! icheckmpc = value of mpc-check parameter option
|
||||||
|
! noline = set to 1 in osolty if no line seacrh should be done in ogetst
|
||||||
|
! icuring = set to 1 if the curing is included for the heat transfer analysis.
|
||||||
|
! ishrink = set to 1 if shrinkage strain is included for mechancial analysis.
|
||||||
|
! ioffsflg = 1 for small displacement beam/shell offsets
|
||||||
|
! = 2 for large displacement beam/shell offsets
|
||||||
|
! isetoff = 0 - do not apply beam/shell offsets
|
||||||
|
! = 1 - apply beam/shell offsets
|
||||||
|
! ioffsetm = min. value of offset flag
|
||||||
|
! inc_incdat = flag to record increment number of a new loadcase in incdat.f
|
||||||
|
! iautspc = flag for AutoSPC option
|
||||||
|
! ibrake = brake squeal in this increment
|
||||||
|
! icbush = set to 1 if cbush elements present in model
|
||||||
|
! istream_input = set to 1 for streaming input calling Marc as library
|
||||||
|
! iprsinp = set to 1 if pressure input, introduced so other variables
|
||||||
|
! such as h could be a function of pressure
|
||||||
|
! ivlsinp = set to 1 if velocity input, introduced so other variables
|
||||||
|
! such as h could be a function of velocity
|
||||||
|
! ipin_m = # of beam element with PIN flag
|
||||||
|
! jgnstr_glb = global control over pre or fast integrated composite shells
|
||||||
|
! imarc_return = Marc return flag for streaming input control
|
||||||
|
! iqvcimp = if non-zero, then the number of QVECT boundary conditions
|
||||||
|
! nqvceid = number of QVECT boundary conditions, where emisivity/absorbtion id entered
|
||||||
|
! istpnx = 1 if to stop at end of increment
|
||||||
|
! imicro1 = 1 if micro1 interface is used
|
||||||
|
! iaxisymm = set to 1 if axisymmetric analysis
|
||||||
|
! jbreakglue = set to 1 if breaking glued option is used
|
||||||
|
! iglstif = 1 if ddm and global stiffness matrix formed (sgi solver 6 or solver9)
|
||||||
|
! jfastasm = 1 do fast assembly using SuperForm code
|
||||||
|
! iwear = set to 1 if wear model, set to 2 if wear model and coordinates updated
|
||||||
|
! iwearcf = set to 1 to store nodal coefficient of friction for wear calculation
|
||||||
|
! imixmeth = set=1 then use nonlinear mixture material - allocate memory
|
||||||
|
! ielcmadyn = flag for magnetodynamics
|
||||||
|
! 0 - electromagnetics using newmark beta
|
||||||
|
! 1 - transient magnetics using backward euler
|
||||||
|
! idinout = flag to control if inside out elements should be deactivated
|
||||||
|
! igena_meth = 0 - generalized alpha parameters depend on whether or not contact
|
||||||
|
! is flagged (dynamic,7)
|
||||||
|
! 10 - generalized alpha parameters are optimized for a contact
|
||||||
|
! analysis (dynamic,8)
|
||||||
|
! 11 - generalized alpha parameters are optimized for an analysis
|
||||||
|
! without contact (dynamic,8)
|
||||||
|
!
|
||||||
|
!***********************************************************************
|
|
@ -0,0 +1,28 @@
|
||||||
|
! reformated to free format
|
||||||
|
!***********************************************************************
|
||||||
|
!
|
||||||
|
! File: creeps.cmn
|
||||||
|
!
|
||||||
|
! MSC.Marc include file
|
||||||
|
!
|
||||||
|
real(pReal) cptim,timinc,timinc_p,timinc_s,timincm,timinc_a,timinc_b,creept
|
||||||
|
integer(pInt) icptim,icfte,icfst,icfeq,icftm,icetem,mcreep,jcreep,icpa,icftmp,icfstr,&
|
||||||
|
icfqcp,icfcpm,icrppr,icrcha,icpb,iicpmt,iicpa
|
||||||
|
real(pReal) time_beg_lcase,time_beg_inc,fractol,time_beg_pst
|
||||||
|
!
|
||||||
|
integer num_creepsr,num_creepsi,num_creeps2r
|
||||||
|
parameter(num_creepsr=40)
|
||||||
|
parameter(num_creepsi=18)
|
||||||
|
parameter(num_creeps2r=4)
|
||||||
|
common/marc_creeps/cptim,timinc,timinc_p,timinc_s,timincm,timinc_a,timinc_b,creept(33),icptim,icfte,icfst,&
|
||||||
|
icfeq,icftm,icetem,mcreep,jcreep,icpa,icftmp,icfstr,icfqcp,icfcpm,icrppr,icrcha,icpb,iicpmt,iicpa
|
||||||
|
common/marc_creeps2/time_beg_lcase,time_beg_inc,fractol,time_beg_pst
|
||||||
|
!
|
||||||
|
! time_beg_lcase time at the beginning of the current load case
|
||||||
|
! time_beg_inc time at the beginning of the current increment
|
||||||
|
! fractol fraction of loadcase or increment time when we
|
||||||
|
! consider it to be finished
|
||||||
|
! time_beg_pst time corresponding to first increment to be
|
||||||
|
! read in from thermal post file for auto step
|
||||||
|
!
|
||||||
|
!***********************************************************************
|
|
@ -0,0 +1,292 @@
|
||||||
|
!********************************************************************
|
||||||
|
! Material subroutine for MSC.Marc Version 0.1
|
||||||
|
!
|
||||||
|
! written by F. Roters, P. Eisenlohr, L. Hantcherli, W.A. Counts
|
||||||
|
! MPI fuer Eisenforschung, Duesseldorf
|
||||||
|
!
|
||||||
|
! last modified: 28.10.2008
|
||||||
|
!********************************************************************
|
||||||
|
! Usage:
|
||||||
|
! - choose material as hypela2
|
||||||
|
! - set statevariable 2 to index of material
|
||||||
|
! - set statevariable 3 to index of texture
|
||||||
|
! - choose output of user variables if desired
|
||||||
|
! - make sure the file "mattex.mpie" exists in the working
|
||||||
|
! directory
|
||||||
|
! - use nonsymmetric option for solver (e.g. direct
|
||||||
|
! profile or multifrontal sparse, the latter seems
|
||||||
|
! to be faster!)
|
||||||
|
! - in case of ddm a symmetric solver has to be used
|
||||||
|
!********************************************************************
|
||||||
|
! Marc subroutines used:
|
||||||
|
! - hypela2
|
||||||
|
! - plotv
|
||||||
|
! - quit
|
||||||
|
!********************************************************************
|
||||||
|
! Marc common blocks included:
|
||||||
|
! - concom: lovl, ncycle, inc, incsub
|
||||||
|
! - creeps: timinc
|
||||||
|
!********************************************************************
|
||||||
|
!
|
||||||
|
include "prec.f90" ! uses nothing else
|
||||||
|
include "debug.f90" ! uses prec
|
||||||
|
include "math.f90" ! uses prec
|
||||||
|
include "IO.f90" ! uses prec, debug, math
|
||||||
|
include "FEsolving.f90" ! uses prec, IO
|
||||||
|
include "mesh.f90" ! uses prec, IO, math, FEsolving
|
||||||
|
include "lattice.f90" ! uses prec, math
|
||||||
|
include "constitutive.f90" ! uses prec, IO, math, lattice, mesh, debug
|
||||||
|
include "crystallite.f90" ! uses prec, debug, constitutive, mesh, math, IO
|
||||||
|
include "CPFEM.f90" ! uses prec, math, mesh, constitutive, FEsolving, debug, lattice, IO, crystallite
|
||||||
|
!
|
||||||
|
!
|
||||||
|
SUBROUTINE hypela2(d,g,e,de,s,t,dt,ngens,n,nn,kcus,matus,ndi,&
|
||||||
|
nshear,disp,dispt,coord,ffn,frotn,strechn,eigvn,ffn1,&
|
||||||
|
frotn1,strechn1,eigvn1,ncrd,itel,ndeg,ndm,&
|
||||||
|
nnode,jtype,lclass,ifr,ifu)
|
||||||
|
!********************************************************************
|
||||||
|
! This is the Marc material routine
|
||||||
|
!********************************************************************
|
||||||
|
!
|
||||||
|
! ************* user subroutine for defining material behavior **************
|
||||||
|
!
|
||||||
|
!
|
||||||
|
! CAUTION : Due to calculation of the Deformation gradients, Stretch Tensors and
|
||||||
|
! Rotation tensors at previous and current states, the analysis can be
|
||||||
|
! computationally expensive. Please use the user subroutine -> hypela
|
||||||
|
! if these kinematic quantities are not needed in the constitutive model
|
||||||
|
!
|
||||||
|
!
|
||||||
|
! IMPORTANT NOTES :
|
||||||
|
!
|
||||||
|
! (1) F,R,U are only available for continuum and membrane elements (not for
|
||||||
|
! shells and beams).
|
||||||
|
!
|
||||||
|
! (2) For total Lagrangian formulation use the -> 'Elasticity,1' card(=
|
||||||
|
! total Lagrange with large disp) in the parameter section of input deck.
|
||||||
|
! For updated Lagrangian formulation use the -> 'Plasticity,3' card(=
|
||||||
|
! update+finite+large disp+constant d) in the parameter section of
|
||||||
|
! input deck.
|
||||||
|
!
|
||||||
|
!
|
||||||
|
! d stress strain law to be formed
|
||||||
|
! g change in stress due to temperature effects
|
||||||
|
! e total elastic strain
|
||||||
|
! de increment of strain
|
||||||
|
! s stress - should be updated by user
|
||||||
|
! t state variables (comes in at t=n, must be updated
|
||||||
|
! to have state variables at t=n+1)
|
||||||
|
! dt increment of state variables
|
||||||
|
! ngens size of stress - strain law
|
||||||
|
! n element number
|
||||||
|
! nn integration point number
|
||||||
|
! kcus(1) layer number
|
||||||
|
! kcus(2) internal layer number
|
||||||
|
! matus(1) user material identification number
|
||||||
|
! matus(2) internal material identification number
|
||||||
|
! ndi number of direct components
|
||||||
|
! nshear number of shear components
|
||||||
|
! disp incremental displacements
|
||||||
|
! dispt displacements at t=n (at assembly, lovl=4) and
|
||||||
|
! displacements at t=n+1 (at stress recovery, lovl=6)
|
||||||
|
! coord coordinates
|
||||||
|
! ncrd number of coordinates
|
||||||
|
! ndeg number of degrees of freedom
|
||||||
|
! itel dimension of F and R, either 2 or 3
|
||||||
|
! nnode number of nodes per element
|
||||||
|
! jtype element type
|
||||||
|
! lclass element class
|
||||||
|
! ifr set to 1 if R has been calculated
|
||||||
|
! ifu set to 1 if strech has been calculated
|
||||||
|
!
|
||||||
|
! at t=n :
|
||||||
|
!
|
||||||
|
! ffn deformation gradient
|
||||||
|
! frotn rotation tensor
|
||||||
|
! strechn square of principal stretch ratios, lambda(i)
|
||||||
|
! eigvn(i,j) i principal direction components for j eigenvalues
|
||||||
|
!
|
||||||
|
! at t=n+1 :
|
||||||
|
!
|
||||||
|
! ffn1 deformation gradient
|
||||||
|
! frotn1 rotation tensor
|
||||||
|
! strechn1 square of principal stretch ratios, lambda(i)
|
||||||
|
! eigvn1(i,j) i principal direction components for j eigenvalues
|
||||||
|
!
|
||||||
|
! The following operation obtains U (stretch tensor) at t=n+1 :
|
||||||
|
!
|
||||||
|
! call scla(un1,0.d0,itel,itel,1)
|
||||||
|
! do 3 k=1,3
|
||||||
|
! do 2 i=1,3
|
||||||
|
! do 1 j=1,3
|
||||||
|
! un1(i,j)=un1(i,j)+dsqrt(strechn1(k))*eigvn1(i,k)*eigvn1(j,k)
|
||||||
|
!1 continue
|
||||||
|
!2 continue
|
||||||
|
!3 continue
|
||||||
|
!
|
||||||
|
use prec, only: pReal,pInt, ijaco
|
||||||
|
use FEsolving
|
||||||
|
use CPFEM, only: CPFEM_general
|
||||||
|
use math, only: invnrmMandel
|
||||||
|
!
|
||||||
|
implicit none
|
||||||
|
!
|
||||||
|
! ** Start of generated type statements **
|
||||||
|
real(pReal) coord, d, de, disp, dispt, dt, e, eigvn, eigvn1, ffn, ffn1
|
||||||
|
real(pReal) frotn, frotn1, g
|
||||||
|
integer(pInt) ifr, ifu, itel, jtype, kcus, lclass, matus, n, ncrd, ndeg
|
||||||
|
integer(pInt) ndi, ndm, ngens, nn, nnode, nshear
|
||||||
|
real(pReal) s, strechn, strechn1, t
|
||||||
|
! ** End of generated type statements **
|
||||||
|
!
|
||||||
|
dimension e(*),de(*),t(*),dt(*),g(*),d(ngens,*),s(*), n(2),coord(ncrd,*),disp(ndeg,*),matus(2),dispt(ndeg,*),ffn(itel,*),&
|
||||||
|
frotn(itel,*),strechn(itel),eigvn(itel,*),ffn1(itel,*),frotn1(itel,*),strechn1(itel),eigvn1(itel,*),kcus(2),&
|
||||||
|
lclass(2)
|
||||||
|
!
|
||||||
|
! Marc common blocks are in fixed format so they have to be reformated to free format (f90)
|
||||||
|
! Beware of changes in newer Marc versions -- these are from 2005r3
|
||||||
|
! concom is needed for inc, subinc, ncycle, lovl
|
||||||
|
include "concom_f90"
|
||||||
|
! creeps is needed for timinc (time increment)
|
||||||
|
include "creeps_f90"
|
||||||
|
!
|
||||||
|
integer(pInt) computationMode,i
|
||||||
|
!
|
||||||
|
if (inc == 0) then
|
||||||
|
cycleCounter = 4
|
||||||
|
else
|
||||||
|
if (theCycle > ncycle .or. theInc /= inc) cycleCounter = 0 ! reset counter for each cutback or new inc
|
||||||
|
if (theCycle /= ncycle .or. theLovl /= lovl) then
|
||||||
|
cycleCounter = cycleCounter+1 ! ping pong
|
||||||
|
outdatedFFN1 = .false.
|
||||||
|
endif
|
||||||
|
endif
|
||||||
|
if (cptim > theTime .or. theInc /= inc) then ! reached convergence
|
||||||
|
lastIncConverged = .true.
|
||||||
|
outdatedByNewInc = .true.
|
||||||
|
endif
|
||||||
|
|
||||||
|
if (mod(cycleCounter,2) /= 0) computationMode = 4 ! recycle
|
||||||
|
if (mod(cycleCounter,4) == 2) computationMode = 3 ! collect
|
||||||
|
if (mod(cycleCounter,4) == 0) computationMode = 2 ! compute
|
||||||
|
if (computationMode == 4 .and. ncycle == 0 .and. .not. lastIncConverged) &
|
||||||
|
computationMode = 6 ! recycle but restore known good consistent tangent
|
||||||
|
if (computationMode == 4 .and. lastIncConverged) then
|
||||||
|
computationMode = 5 ! recycle and record former consistent tangent
|
||||||
|
lastIncConverged = .false.
|
||||||
|
endif
|
||||||
|
if (computationMode == 2 .and. outdatedByNewInc) then
|
||||||
|
computationMode = 1 ! compute and age former results
|
||||||
|
outdatedByNewInc = .false.
|
||||||
|
endif
|
||||||
|
if (computationMode == 2 .and. outdatedFFN1) then
|
||||||
|
computationMode = 4 ! return odd results to force new vyvle
|
||||||
|
endif
|
||||||
|
|
||||||
|
|
||||||
|
theTime = cptim ! record current starting time
|
||||||
|
theInc = inc ! record current increment number
|
||||||
|
theCycle = ncycle ! record current cycle count
|
||||||
|
theLovl = lovl ! record current lovl
|
||||||
|
|
||||||
|
call CPFEM_general(computationMode,ffn,ffn1,t(1),timinc,n(1),nn,s,mod(cycleCounter-4,4_pInt*ijaco)==0,d,ngens)
|
||||||
|
|
||||||
|
|
||||||
|
! Mandel: 11, 22, 33, SQRT(2)*12, SQRT(2)*23, SQRT(2)*13
|
||||||
|
! Marc: 11, 22, 33, 12, 23, 13
|
||||||
|
forall(i=1:ngens) d(1:ngens,i) = invnrmMandel(i)*d(1:ngens,i)*invnrmMandel(1:ngens)
|
||||||
|
s(1:ngens) = s(1:ngens)*invnrmMandel(1:ngens)
|
||||||
|
if(symmetricSolver) d(1:ngens,1:ngens) = 0.5_pReal*(d(1:ngens,1:ngens)+transpose(d(1:ngens,1:ngens)))
|
||||||
|
return
|
||||||
|
|
||||||
|
END SUBROUTINE
|
||||||
|
!
|
||||||
|
|
||||||
|
SUBROUTINE plotv(v,s,sp,etot,eplas,ecreep,t,m,nn,layer,ndi,nshear,jpltcd)
|
||||||
|
!********************************************************************
|
||||||
|
! This routine sets user defined output variables for Marc
|
||||||
|
!********************************************************************
|
||||||
|
!
|
||||||
|
! select a variable contour plotting (user subroutine).
|
||||||
|
!
|
||||||
|
! v variable
|
||||||
|
! s (idss) stress array
|
||||||
|
! sp stresses in preferred direction
|
||||||
|
! etot total strain (generalized)
|
||||||
|
! eplas total plastic strain
|
||||||
|
! ecreep total creep strain
|
||||||
|
! t current temperature
|
||||||
|
! m element number
|
||||||
|
! nn integration point number
|
||||||
|
! layer layer number
|
||||||
|
! ndi (3) number of direct stress components
|
||||||
|
! nshear (3) number of shear stress components
|
||||||
|
!
|
||||||
|
!********************************************************************
|
||||||
|
use prec, only: pReal,pInt
|
||||||
|
use CPFEM, only: CPFEM_results, CPFEM_Nresults
|
||||||
|
use constitutive, only: constitutive_maxNresults
|
||||||
|
use mesh, only: mesh_FEasCP
|
||||||
|
implicit none
|
||||||
|
!
|
||||||
|
real(pReal) s(*),etot(*),eplas(*),ecreep(*),sp(*)
|
||||||
|
real(pReal) v, t(*)
|
||||||
|
integer(pInt) m, nn, layer, ndi, nshear, jpltcd
|
||||||
|
!
|
||||||
|
! assign result variable
|
||||||
|
v=CPFEM_results(mod(jpltcd-1_pInt, CPFEM_Nresults+constitutive_maxNresults)+1_pInt,&
|
||||||
|
(jpltcd-1_pInt)/(CPFEM_Nresults+constitutive_maxNresults)+1_pInt,&
|
||||||
|
nn, mesh_FEasCP('elem', m))
|
||||||
|
return
|
||||||
|
END SUBROUTINE
|
||||||
|
!
|
||||||
|
!
|
||||||
|
! subroutine utimestep(timestep,timestepold,icall,time,timeloadcase)
|
||||||
|
!********************************************************************
|
||||||
|
! This routine modifies the addaptive time step of Marc
|
||||||
|
!********************************************************************
|
||||||
|
! use prec, only: pReal,pInt
|
||||||
|
! use CPFEM, only : CPFEM_timefactor_max
|
||||||
|
! implicit none
|
||||||
|
!
|
||||||
|
! real(pReal) timestep, timestepold, time,timeloadcase
|
||||||
|
! integer(pInt) icall
|
||||||
|
!
|
||||||
|
! user subroutine for modifying the time step in auto step
|
||||||
|
!
|
||||||
|
! timestep : the current time step as suggested by marc
|
||||||
|
! to be modified in this routine
|
||||||
|
! timestepold : the current time step before it was modified by marc
|
||||||
|
! icall : =1 for setting the initial time step
|
||||||
|
! =2 if this routine is called during an increment
|
||||||
|
! =3 if this routine is called at the beginning
|
||||||
|
! of the increment
|
||||||
|
! time : time at the start of the current increment
|
||||||
|
! timeloadcase: time period of the current load case
|
||||||
|
!
|
||||||
|
! it is in general not recommended to increase the time step
|
||||||
|
! during the increment.
|
||||||
|
! this routine is called right after the time step has (possibly)
|
||||||
|
! been updated by marc.
|
||||||
|
!
|
||||||
|
! user coding
|
||||||
|
! reduce timestep during increment in case mpie_timefactor is too large
|
||||||
|
! if(icall==2_pInt) then
|
||||||
|
! if(mpie_timefactor_max>1.25_pReal) then
|
||||||
|
! timestep=min(timestep,timestepold*0.8_pReal)
|
||||||
|
! end if
|
||||||
|
! return
|
||||||
|
! modify timestep at beginning of new increment
|
||||||
|
! else if(icall==3_pInt) then
|
||||||
|
! if(mpie_timefactor_max<=0.8_pReal) then
|
||||||
|
! timestep=min(timestep,timestepold*1.25_pReal)
|
||||||
|
! else if (mpie_timefactor_max<=1.0_pReal) then
|
||||||
|
! timestep=min(timestep,timestepold/mpie_timefactor_max)
|
||||||
|
! else if (mpie_timefactor_max<=1.25_pReal) then
|
||||||
|
! timestep=min(timestep,timestepold*1.01_pReal)
|
||||||
|
! else
|
||||||
|
! timestep=min(timestep,timestepold*0.8_pReal)
|
||||||
|
! end if
|
||||||
|
! end if
|
||||||
|
! return
|
||||||
|
! end
|
Loading…
Reference in New Issue