10   private :: update_abort
 
   30     type (hecmwST_matrix)       :: hecMAT
 
   31     type (hecmwST_local_mesh)   :: hecMESH
 
   32     type (fstr_solid)           :: fstrSOLID
 
   33     real(kind=kreal),
intent(in) :: time      
 
   34     real(kind=kreal),
intent(in) :: tincr     
 
   35     integer, 
intent(in)         :: iter
 
   37     integer(kind=kint) :: nodLOCAL(fstrSOLID%max_ncon)
 
   38     real(kind=kreal)   :: ecoord(3, fstrsolid%max_ncon)
 
   39     real(kind=kreal)   :: thick, thick0(6)
 
   40     integer(kind=kint) :: ndof, itype, is, iE, ic_type, nn, icel, iiS, i, j
 
   42     real(kind=kreal)   :: total_disp(6, fstrsolid%max_ncon), du(6, fstrsolid%max_ncon), ddu(6, fstrsolid%max_ncon)
 
   43     real(kind=kreal)   :: tt(fstrsolid%max_ncon), tt0(fstrsolid%max_ncon), ttn(fstrsolid%max_ncon)
 
   44     real(kind=kreal)   :: qf(fstrsolid%max_ncon*6), coords(3, 3)
 
   45     integer            :: isect, ihead, cdsys_ID
 
   46     integer            :: ndim, initt
 
   48     real(kind=kreal), 
optional :: strainenergy
 
   49     real(kind=kreal) :: tmp
 
   50     real(kind=kreal)   :: ddaux(3,3)
 
   53     fstrsolid%QFORCE=0.0d0
 
   80     do itype = 1, hecmesh%n_elem_type
 
   81       is = hecmesh%elem_type_index(itype-1)+1
 
   82       ie = hecmesh%elem_type_index(itype  )
 
   83       ic_type= hecmesh%elem_type_item(itype)
 
   84       if (hecmw_is_etype_link(ic_type)) cycle
 
   85       if (hecmw_is_etype_patch(ic_type)) cycle
 
   98         iis = hecmesh%elem_node_index(icel-1)
 
   99         nn = hecmesh%elem_node_index(icel)-iis
 
  103         do j = 1, 
size(fstrsolid%elements(icel)%gausses)
 
  104           thick0(1:6) = fstrsolid%elements(icel)%gausses(j)%stress(1:6)
 
  105           thick = thick + dsqrt(dot_product(thick0(1:6),thick0(1:6)))
 
  107         if( thick < 1.d-10 ) 
then 
  108           do j = 1, 
size(fstrsolid%elements(icel)%gausses)
 
  109             if( 
associated(fstrsolid%elements(icel)%gausses(j)%fstatus) ) &
 
  110               &  fstrsolid%elements(icel)%gausses(j)%fstatus = 0.d0
 
  115           nodlocal(j) = hecmesh%elem_node_item (iis+j)
 
  117             ecoord(i,j) = hecmesh%node(3*nodlocal(j)+i-3)
 
  120             ddu(i,j) = hecmat%X(ndof*nodlocal(j)+i-ndof)
 
  121             du(i,j)  = fstrsolid%dunode(ndof*nodlocal(j)+i-ndof)
 
  122             total_disp(i,j) = fstrsolid%unode(ndof*nodlocal(j)+i-ndof)
 
  125           if( fstrsolid%TEMP_ngrp_tot > 0 .or. fstrsolid%TEMP_irres > 0 ) 
then 
  126             if( iselastoplastic(fstrsolid%elements(icel)%gausses(1)%pMaterial%mtype) .or. &
 
  127                 fstrsolid%elements(icel)%gausses(1)%pMaterial%mtype == norton ) 
then 
  128               tt0(j)=fstrsolid%last_temp( nodlocal(j) )
 
  131               if( hecmesh%hecmw_flag_initcon == 1 ) tt0(j) = hecmesh%node_init_val_item(nodlocal(j))
 
  132               if( initt>0 ) tt0(j) = 
g_initialcnd(initt)%realval(nodlocal(j))
 
  134             ttn(j) = fstrsolid%last_temp( nodlocal(j) )
 
  135             tt(j)  = fstrsolid%temperature( nodlocal(j) )
 
  139         isect = hecmesh%section_ID(icel)
 
  140         ihead = hecmesh%section%sect_R_index(isect-1)
 
  141         thick = hecmesh%section%sect_R_item(ihead+1)
 
  142         cdsys_id = hecmesh%section%sect_orien_ID(isect)
 
  143         if( cdsys_id > 0 ) 
call get_coordsys(cdsys_id, hecmesh, fstrsolid, coords)
 
  146         if( getspacedimension( ic_type ) == 2 ) thick = 1.0d0
 
  148         if( ic_type == 241 .or. ic_type == 242 .or. ic_type == 231 .or. ic_type == 232 .or. ic_type == 2322 ) 
then 
  149           call update_c2( ic_type,nn,ecoord(1:3,1:nn),fstrsolid%elements(icel)%gausses(:), &
 
  150             thick,fstrsolid%elements(icel)%iset,                             &
 
  151             total_disp(1:2,1:nn), ddu(1:2,1:nn), qf(1:nn*ndof),              &
 
  152             tt(1:nn), tt0(1:nn), ttn(1:nn)  )
 
  154         else if( ic_type == 301 ) 
then 
  155           call update_c1( ic_type,nn,ecoord(:,1:nn), thick, total_disp(1:3,1:nn), du(1:3,1:nn), &
 
  156             qf(1:nn*ndof),fstrsolid%elements(icel)%gausses(:) )
 
  158         else if( ic_type == 361 ) 
then 
  159           if( fstrsolid%sections(isect)%elemopt361 == 
kel361fi ) 
then  
  160             call update_c3( ic_type, nn, ecoord(:,1:nn), total_disp(1:3,1:nn), du(1:3,1:nn), cdsys_id, coords, &
 
  161               qf(1:nn*ndof), fstrsolid%elements(icel)%gausses(:), iter, time, tincr, tt(1:nn), tt0(1:nn), ttn(1:nn)  )
 
  162           else if( fstrsolid%sections(isect)%elemopt361 == 
kel361bbar ) 
then  
  163             call update_c3d8bbar( ic_type, nn, ecoord(:,1:nn), total_disp(1:3,1:nn), du(1:3,1:nn), cdsys_id, coords,    &
 
  164               qf(1:nn*ndof), fstrsolid%elements(icel)%gausses(:), iter, time, tincr, tt(1:nn), tt0(1:nn), ttn(1:nn)  )
 
  165           else if( fstrsolid%sections(isect)%elemopt361 == 
kel361ic ) 
then  
  166             call update_c3d8ic( ic_type,nn,ecoord(:,1:nn), total_disp(1:3,1:nn), du(1:3,1:nn), ddu(1:3,1:nn), cdsys_id, coords,&
 
  167               qf(1:nn*ndof), fstrsolid%elements(icel)%gausses(:), iter, time, tincr, &
 
  168               fstrsolid%elements(icel)%aux, ddaux(1:3,1:3), tt(1:nn), tt0(1:nn), ttn(1:nn) )
 
  169             fstrsolid%elements(icel)%aux(1:3,1:3) = fstrsolid%elements(icel)%aux(1:3,1:3) + ddaux(1:3,1:3)
 
  170           else if( fstrsolid%sections(isect)%elemopt361 == 
kel361fbar ) 
then  
  171             call update_c3d8fbar( ic_type, nn, ecoord(:,1:nn), total_disp(1:3,1:nn), du(1:3,1:nn), cdsys_id, coords,    &
 
  172               qf(1:nn*ndof), fstrsolid%elements(icel)%gausses(:), iter, time, tincr, tt(1:nn), tt0(1:nn), ttn(1:nn)  )
 
  175         else if (ic_type == 341 .or. ic_type == 351 .or. ic_type == 342 .or. ic_type == 352 .or. ic_type == 362 ) 
then 
  176           if( ic_type==341 .and. fstrsolid%sections(isect)%elemopt341 == 
kel341sesns ) cycle 
 
  177           call update_c3( ic_type, nn, ecoord(:,1:nn), total_disp(1:3,1:nn), du(1:3,1:nn), cdsys_id, coords, &
 
  178             qf(1:nn*ndof), fstrsolid%elements(icel)%gausses(:), iter, time, tincr, tt(1:nn), tt0(1:nn), ttn(1:nn)  )
 
  180         else if( ic_type == 511) 
then 
  181           call update_connector( ic_type,nn,ecoord(:,1:nn), total_disp(1:3,1:nn), du(1:3,1:nn), &
 
  182             qf(1:nn*ndof),fstrsolid%elements(icel)%gausses(:), tincr )
 
  184         else if( ic_type == 611) 
then 
  185           if( 
fstrpr%nlgeom ) 
call update_abort( ic_type, 2 )
 
  186           CALL updatest_beam(ic_type, nn, ecoord, total_disp(1:6,1:nn), du(1:6,1:nn), &
 
  187                  &   hecmesh%section%sect_R_item(ihead+1:), fstrsolid%elements(icel)%gausses(:), qf(1:nn*ndof))
 
  189         else if( ic_type == 641 ) 
then 
  190           if( 
fstrpr%nlgeom ) 
call update_abort( ic_type, 2 )
 
  191           call updatest_beam_641(ic_type, nn, ecoord, total_disp(1:ndof,1:nn), du(1:ndof,1:nn), &
 
  192             &    fstrsolid%elements(icel)%gausses(:), hecmesh%section%sect_R_item(ihead+1:), qf(1:nn*ndof))
 
  194         else if( ( ic_type == 741 ) .or. ( ic_type == 743 ) .or. ( ic_type == 731 ) ) 
then 
  195           if( 
fstrpr%nlgeom ) 
call update_abort( ic_type, 2 )
 
  196           call updatest_shell_mitc(ic_type, nn, ndof, ecoord(1:3, 1:nn), total_disp(1:ndof,1:nn), du(1:ndof,1:nn), &
 
  197             &              fstrsolid%elements(icel)%gausses(:), qf(1:nn*ndof), thick, 0)
 
  199         else if( ic_type == 761 ) 
then    
  200           if( 
fstrpr%nlgeom ) 
call update_abort( ic_type, 2 )
 
  201           call updatest_shell_mitc33(731, 3, 6, ecoord(1:3, 1:3), total_disp(1:ndof,1:nn), du(1:ndof,1:nn), &
 
  202             &              fstrsolid%elements(icel)%gausses(:), qf(1:nn*ndof), thick, 2)
 
  204         else if( ic_type == 781 ) 
then    
  205           if( 
fstrpr%nlgeom ) 
call update_abort( ic_type, 2 )
 
  206           call updatest_shell_mitc33(741, 4, 6, ecoord(1:3, 1:4), total_disp(1:ndof,1:nn), du(1:ndof,1:nn), &
 
  207             &              fstrsolid%elements(icel)%gausses(:), qf(1:nn*ndof), thick, 1)
 
  209         else if ( ic_type == 3414 ) 
then 
  210           if(fstrsolid%elements(icel)%gausses(1)%pMaterial%mtype /= incomp_newtonian) &
 
  211             &  
call update_abort( ic_type, 3, fstrsolid%elements(icel)%gausses(1)%pMaterial%mtype )
 
  213             ( ic_type, nn, ecoord(:,1:nn), total_disp(1:4,1:nn), du(1:4,1:nn), &
 
  214             fstrsolid%elements(icel)%gausses(:) )
 
  217         else if ( ic_type == 881 .or. ic_type == 891 ) 
then   
  218           call update_c3_sesns( ic_type, nn, nodlocal, ecoord(:,1:nn), total_disp(1:3,1:nn), du(1:3,1:nn), cdsys_id, coords, &
 
  219             qf(1:nn*ndof), fstrsolid%elements(icel)%gausses(:), time, tincr, tt(1:nn), tt0(1:nn), ttn(1:nn)  )
 
  222           write(*, *) 
'###ERROR### : Element type not supported for nonlinear static analysis' 
  223           write(*, *) 
' ic_type = ', ic_type
 
  224           call hecmw_abort(hecmw_comm_get_comm())
 
  230           call update_dummy( ndof, nn, ecoord(:,1:nn), total_disp(1:3,1:nn), &
 
  231             &  du(1:3,1:nn), qf(1:nn*ndof), fstrsolid%elements(icel) )
 
  239             fstrsolid%QFORCE(ndof*(nodlocal(j)-1)+i) = fstrsolid%QFORCE(ndof*(nodlocal(j)-1)+i)+qf(ndof*(j-1)+i)
 
  244         if(
present(strainenergy))
then 
  245           ndim = getspacedimension( fstrsolid%elements(icel)%etype )
 
  248               tmp = 0.5d0*( fstrsolid%elements(icel)%equiForces(ndim*(j-1)+i)+qf(ndim*(j-1)+i) )*ddu(i,j)
 
  250               strainenergy = strainenergy+tmp
 
  251               fstrsolid%elements(icel)%equiForces(ndim*(j-1)+i) = qf(ndim*(j-1)+i)
 
  264     call hecmw_update_r(hecmesh,fstrsolid%QFORCE,hecmesh%n_node, ndof)
 
  275     type(hecmwst_local_mesh) :: hecmesh
 
  277     real(kind=kreal) :: tincr
 
  278     integer(kind=kint) :: itype, is, iE, ic_type, icel, ngauss, i
 
  280     if( 
associated( fstrsolid%temperature ) ) 
then 
  281       do i = 1, hecmesh%n_node
 
  282         fstrsolid%last_temp(i) = fstrsolid%temperature(i)
 
  286     do itype = 1, hecmesh%n_elem_type
 
  287       is = hecmesh%elem_type_index(itype-1) + 1
 
  288       ie = hecmesh%elem_type_index(itype  )
 
  289       ic_type= hecmesh%elem_type_item(itype)
 
  290       if( ic_type == 301 ) ic_type = 111
 
  291       if( hecmw_is_etype_link(ic_type) ) cycle
 
  292       if( hecmw_is_etype_patch(ic_type) ) cycle
 
  294       ngauss = numofquadpoints( ic_type )
 
  296         if( iselastoplastic( fstrsolid%elements(icel)%gausses(1)%pMaterial%mtype ) ) 
then 
  300         elseif( fstrsolid%elements(icel)%gausses(1)%pMaterial%mtype == norton ) 
then 
  301           if( tincr>0.d0 ) 
then 
  306         elseif( isviscoelastic( fstrsolid%elements(icel)%gausses(1)%pMaterial%mtype ) ) 
then 
  307           if( tincr > 0.d0 ) 
then 
  315           fstrsolid%elements(icel)%gausses(i)%strain_bak = fstrsolid%elements(icel)%gausses(i)%strain
 
  316           fstrsolid%elements(icel)%gausses(i)%stress_bak = fstrsolid%elements(icel)%gausses(i)%stress
 
  317           fstrsolid%elements(icel)%gausses(i)%strain_energy_bak = fstrsolid%elements(icel)%gausses(i)%strain_energy
 
  322     do i = 1, hecmesh%n_node
 
  323       fstrsolid%QFORCE_bak(i) = fstrsolid%QFORCE(i)
 
  328   subroutine update_abort( ic_type, flag, mtype )
 
  329     integer(kind=kint), 
intent(in) :: ic_type
 
  330     integer(kind=kint), 
intent(in) :: flag
 
  331     integer(kind=kint), 
intent(in), 
optional :: mtype
 
  334       write(*,*) 
'###ERROR### : Element type not supported for static analysis' 
  335     else if( flag == 2 ) 
then 
  336       write(*,*) 
'###ERROR### : Element type not supported for nonlinear static analysis' 
  337     else if( flag == 3 ) 
then 
  338       write(*,*) 
'###ERROR### : This element is not supported for this material' 
  340     write(*,*) 
' ic_type = ', ic_type
 
  341     if( 
present(mtype) ) 
write(*,*) 
' mtype = ', mtype
 
  342     call hecmw_abort(hecmw_comm_get_comm())
 
This module provide functions for elastoplastic calculation.
 
subroutine, public updateepstate(gauss)
Clear elatoplastic state.
 
This module defined elemact data and function.
 
subroutine update_dummy(ndof, nn, ecoord, u, du, qf, element)
 
integer, parameter kelact_inactive
 
This module provides function to calculate to do updates.
 
subroutine fstr_updatestate(hecMESH, fstrSOLID, tincr)
Update elastiplastic status.
 
subroutine fstr_updatenewton(hecMESH, hecMAT, fstrSOLID, time, tincr, iter, strainEnergy)
Update displacement, stress, strain and internal forces.
 
This module defines common data and basic structures for analysis.
 
integer(kind=kint), parameter kel361bbar
 
subroutine get_coordsys(cdsys_ID, hecMESH, fstrSOLID, coords)
This subroutine fetch coords defined by local coordinate system.
 
integer(kind=kint), parameter kel341sesns
 
integer(kind=kint), parameter kel361fi
 
integer(kind=kint), parameter kel361ic
 
type(fstr_param), target fstrpr
GLOBAL VARIABLE INITIALIZED IN FSTR_SETUP.
 
integer(kind=kint), parameter kel361fbar
 
type(tinitialcondition), dimension(:), pointer, save g_initialcnd
 
This modules just summarizes all modules used in static analysis.
 
This module provides functions for creep calculation.
 
subroutine updateviscostate(gauss)
Update viscoplastic state.
 
This module provides functions for viscoelastic calculation.
 
subroutine updateviscoelasticstate(gauss)
Update viscoplastic state.