[Paraview] paraview file-formats

SamuelKey samuelkey at comcast.net
Wed Aug 30 21:47:01 EDT 2006


Greetings Adrian,

I suggest you consider using EnSight-formatted binary files. See 
attached example of a Fortran-90 language example that writes such data. 
A complete description of EnSight formats can be found in Chapter 11 of 
the EnSight Online User Manual.
https://www.ensight.com/downloads/doc_download-47.html

Each cell variable and each point variable is written to its own file.  
If you learn how to use this format's "part" concept
you can subdivide your mesh into any number of individually displayable 
pieces based on any criterion you chose.  The Extract Datasets filter 
will then let select a subset of 1, 2 or many parts for 
closer/individual examination.

Regards,

Sam Key


Adrian Magda wrote:
> Hi all,
>
> I use paraview to postprocess the results obtained with our home-made CFD
> software using the VTK file format.
>
> Due to the size of the files (even if i write binary) i need to split the
> results in several files.
>
> 1 - It is necessary to write the geometry in each VTK file (for 
> example one
> contains the temp, other the vel & turbulence, other the species and 
> so on)
> ? If i write a  different file for each property or group of flow
> properties  it seems i can not  color a isosurface of a flow variable
> according to another flow properties from another file (i.e. temp 
> isosurface
> with co2 concentrations)
>
> 2 - What file-type i need to use in order to write the geometry only once
> and have the flow properties in different vtk files?
>
> Thanks a lot
>
> Adi
>
> ------------------------------------------------------------------------
>
> _______________________________________________
> ParaView mailing list
> ParaView at paraview.org
> http://www.paraview.org/mailman/listinfo/paraview
>   
-------------- next part --------------
      SUBROUTINE WRITE_TO_ESG_RESULTS_FILE
!!
!! Copyright (c) by FMA Development, LLC, 28-SEP-2005
!!
!! Purpose: Write results to EnSight Gold C-binary files suitable for
!! EnSight or ParaView postprocessing. Note: EnSight's "UNFORMATTED"
!! Fortran file option (Fortran Binary) could not be used because 
!! ParaView does not support it. See the "FORM = 'BINARY'" directive
!! in the OPEN statements contained in this file.
!!
      USE shared_common_data
!!
!! The complete simulation data set.
!!
      USE indx_;           USE node_;           USE input_function_;
      USE beam_;           USE coord_;          USE sliding_interface_;
      USE value_;          USE force_;          USE nodal_constraints_;
      USE hexah_;          USE penta_;          USE nonreflecting_bc_;
      USE tetra_;          USE lsold_;          USE nodal_point_mass_;
      USE membq_;          USE membt_;          USE rigid_body_mass_;
      USE truss_;          USE platq_;          USE state_variables_;
      USE platt_;          USE motion_;         USE enumerated_sets_;
      USE spring_;         USE damper_;         USE displacement_bc_;
      USE stress_;         USE segment_;        USE contact_surface_;
      USE tied_bc_;        USE results_;        USE relink_scratch_;
      USE gauge1d_;        USE gauge2d_;        USE rigid_wall_bc_;
      USE gauge3d_;        USE massprop_;       USE include_file_;
      USE material_;       USE layering_;       USE sliding_node_;
      USE force_bc_;       USE node_set_;       USE contact_node_;
      USE nrbc_data_;      USE spring_bc_;      USE periodic_bc_;
      USE damper_bc_;      USE spot_weld_;      USE pressure_bc_;
      USE qa_record_;      USE plate_pair_;     USE segment_set_;
      USE body_force_;     USE section_2d_;     USE element_set_;
      USE section_1d_;     USE rigid_body_;     USE velocity_ic_;
      USE section_3d_;     USE extreme_value_;  USE centrifugal_force_;
      USE location_;       USE mean_stress_;    USE output_;         
      USE wedge_;          USE energy_flow_;    USE pyramid_;

      USE precision_

      IMPLICIT    REAL(RTYPE) ( A-H, O-Z )
      IMPLICIT INTEGER(ITYPE) ( I-N      )
!!
!! Local variables.
      INTEGER(ITYPE), SAVE :: TStep = 0   ! Local count of time steps written
      INTEGER(ITYPE), SAVE :: MStep       ! Max number time steps requested/written
      CHARACTER(2)         :: LEGO        ! Used in Pathf95-unique procedure ASSIGN
      CHARACTER(8)         :: MEGO        ! Buffer for Material MatID output
      CHARACTER(5)         :: NEGO        ! Buffer for number-of-time-steps
      CHARACTER(80)        :: CBUFFER     ! EnSight Gold character const buffer

      REAL(RTYPE), DIMENSION(:), ALLOCATABLE, SAVE :: TIME_VALUES

      INTEGER              :: NERROR
      INTEGER(ITYPE), SAVE :: NUMXL
      INTEGER(ITYPE)       :: NTUPLE(8)
      INTEGER(ITYPE)       :: WedgeID
      LOGICAL              :: IOERROR
      LOGICAL,    EXTERNAL :: NEXT_NP_ID
      LOGICAL,    EXTERNAL :: NEXT_SEG_ID
      LOGICAL,        SAVE :: FIRST = .TRUE.
      LOGICAL,        SAVE :: INSERT_C_BINARY_QUE = .TRUE. 

      INTEGER(ITYPE), SAVE :: MHX,MPX,MPY,MTX,MM3,MP3,MM4,MP4,MTR
!!
!! The following arrays are used to extract esg/vtk "parts" based on
!! material models.  For material "M," I = ELEMENTS_AND_NODES_USED(M) 
!! fills the following arrays with indecies.
!!
      INTEGER(ITYPE), DIMENSION(:), ALLOCATABLE, SAVE :: NELUSED ! Elems used by material M
      INTEGER(ITYPE), DIMENSION(:), ALLOCATABLE, SAVE :: NPTUSED ! Nodes used by material M
      INTEGER(ITYPE), DIMENSION(:), ALLOCATABLE, SAVE :: NPTNOWI ! Node's ESG output-index
!!
!! Contained functions
!!    INTEGER :: ELEMENTS_AND_NODES_USED ! Gens material access arrays
!!    CHAR(8) :: ELEMENT_ESG_TYPE        ! Rtns EnSight Gold elem type
!!    INTEGER :: ELEMENT_N_TUPLE         ! Rtns connectivity count
!!    INTEGER :: MatID_DATA              ! Rtns material ID MatID
!!    INTEGER :: EleID_DATA              ! Rtns element ID  EleID
!!    REAL    :: STRESS_DATA             ! Rtns element stress
!!    REAL    :: BULK_STRAIN             ! Rtns element bulk strain
!!    REAL    :: STRAIN_ENERGY_DENSITY   ! Rtns element int energy
!!    REAL    :: PRESSURE                ! Rtns element pressure
!!    REAL    :: EFFECTIVE_STRESS        ! Rtns element eff stress
!!    INTEGER :: SEGMENTS_AND_NODES_USED ! Gens segment access arrays
!!    INTEGER :: SEGMENT_N_TUPLE         ! Rtns connectivity count
!!    INTEGER :: WEDGES_AND_NODES_USED   ! Gens wedge access arrays
!!    INTEGER :: WEDGE_N_TUPLE           ! Rtns connectivity count
!!    REAL    :: WEDGE_NP_COORDINATE     ! Rtns wedge np coordinates
!!    REAL    :: WEDGE_NP_VELOCITY       ! Rtns wedge np velocities
!!
!!############################################################################
!! 1. CONSTANT FOR ALL-TIME DATA
!! "Create one-time geometry for subsequent ESG Data Files to reference."
!! This file will contain the undeformed configuration. It is generated the
!! first time this routine is called and contains the current velocities. 
!!
!! =========FIRST=============================================================
!! This file will contain the the undeformed mesh colored by material and 
!! have current velocity at the nodal points for examination.
!!
      IF (FIRST) THEN

#ifdef LANGUAGE_FORTRAN90
!!
!! ASSIGN is a PathScale pathf95 compiler unique procedure.
!! It causes the unformatted output to be Big_Endian (mips).
!! 
        CALL ASSIGN( "assign -N mips p:fmaego%", NERROR )
        IF (NERROR .NE. 0) THEN
          WRITE (MSG1,'(I8)') NERROR
          CALL USER_MESSAGE
     &      (
     &      MSGL//"FATAL"//
     &      MSGL//"WRITE_TO_ESG_RESULTS_FILE.999.99"//
     &      MSGL//"Call To pathf95 ASSIGN Failed."//
     &      MSGL//"Call PathScale. Error Number:"//MSG1
     &      )
        ENDIF
#endif
!!
!! Allocate storage for the number of time steps expected from the ESGFILE
!! input record entries. If the requested time increment is zero, limit the
!! number of time steps.
!!
      Tbgn = ESG_RESULTS_FILE%Begin
      Tinc = ESG_RESULTS_FILE%Delta
      Tend = MIN( TIME%Stop, ESG_RESULTS_FILE%End )

      IF (Tinc .GT. ZERO) THEN
        MStep = NINT( (Tend-Tbgn)/Tinc ) + 1
      ELSE
        MStep = 101
        WRITE (MSG1,'(I8)') MStep
        CALL USER_MESSAGE
     &    (
     &    MSGL//'INFORM'//
     &    MSGL//'WRITE_ESG_DATA_FILE.010.01'//
     &    MSGL//'Since "ESGFILE DELTA=0.0" Was Input And'//
     &    MSGL//'Due To The Nature Of The EnSight Gold Case'//
     &    MSGL//'File Format, The Number Of Time Steps Was'//
     &    MSGL//'Capped At'//MSG1//' Steps. Try A Non-Zero'//
     &    MSGL//'DELTA Value To Allow A Non-Zero Divide.'
     &    )
      ENDIF

      ALLOCATE ( TIME_VALUES(1:MStep) )
!!
!! Total finite element count. Note: Only one qudrilateral is produced for
!! each membrane (NUMM3, NUMM4) and shell (NUMP3, NUMP4) finite element.
!!
        NUMXL = NUMHX + NUMPX + NUMPY + NUMTX + NUMM3 + NUMP3 + NUMM4 + NUMP4 + NUMTR
!!
!! Count number of elements using each material model. Note: There may be more
!! material models specified then are actually referenced by elements in the
!! mesh. Hence, the effort to isolate materials that have a null usage count.
!!
        MXSNW = MAX ( NUMXL, NUMSG, NUMNP, NUMWX )
        ALLOCATE ( NELUSED(1:MXSNW), NPTUSED(1:NUMNP), NPTNOWI(1:NUMNP) )
!!
!! Define element and node counts for MATERIAL(*)%NElems and MATERIAL(*)%NNodes
!!
        DO M = 1,NUMMT
          I = ELEMENTS_AND_NODES_USED(M)
        ENDDO
!!
!! Open one-time geometry file fmaego.mesh.geom (EnSight Gold output).
!!
        IOERROR = .TRUE.
        OPEN
     &    (
     &    UNIT    =  IO_UNIT%LEGO,
     &    FILE    = 'fmaego.mesh.geom',
     &    STATUS  = 'UNKNOWN',
#ifdef _G95_
     &    FORM    = 'UNFORMATTED', ACCESS='STREAM',  !  G95
#endif
#ifdef _CVF_NT_
     &    FORM    = 'BINARY', CONVERT='BIG_ENDIAN',  !  CVF
#endif
#ifdef LANGUAGE_FORTRAN90
     &    FORM    = 'BINARY',                        !  Pathf95
#endif
     &    ERR     =  100
     &    )
        IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 100    IF (IOERROR) THEN
          CALL USER_MESSAGE
     &      (
     &      MSGL//'FATAL'//
     &      MSGL//'WRITE_ESG_DATA_FILE.001.00'//
     &      MSGL//'Unable To Execute OPEN On: '//'fmaego.mesh.geom'
     &      )
        ELSE
!!
!! Initialize sequential Material part counter.
!!
          MPart = 0
!!
!! Initialize EnSight Gold C-binary data file (see FORM='UNFORMATTED' above).
!!
          CBUFFER = "C Binary"                         ! EnSight reader format que.
          WRITE (IO_UNIT%LEGO) CBUFFER
          CBUFFER = TRIM(JOB_ID_RECORD%CURRENT%TITLE)  ! User's job title record.
          WRITE (IO_UNIT%LEGO) CBUFFER
          CBUFFER = " "                                ! (unused subtitle)
          WRITE (IO_UNIT%LEGO) CBUFFER
          CBUFFER = "node id given"                    ! Expect to read nodal ID's
          WRITE (IO_UNIT%LEGO) CBUFFER
          CBUFFER = "element id given"                 ! Expect to read element ID's
          WRITE (IO_UNIT%LEGO) CBUFFER
!!
!! Turn individual material domains into esg/vtk "parts."
!!
          DO M = 1,NUMMT
!!
!! The function ELEMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! elements-used array NELUSED and vtk output-index translation array NPTNOWI for 
!! this material. If material M is used, the function returns 1, otherwise 0. 
!!
            IF (ELEMENTS_AND_NODES_USED(M) .GT. 0) THEN

              WRITE (MEGO,'(I8.8)') MATERIAL(M)%MatID

              NNodes = MATERIAL(M)%NNodes
              NElems = MATERIAL(M)%NElems

              CBUFFER = "part"
              WRITE (IO_UNIT%LEGO) CBUFFER
              MPart = Mpart + 1
              WRITE (IO_UNIT%LEGO) MPart
              CBUFFER = "Part with MatID: "//MEGO
              WRITE (IO_UNIT%LEGO) CBUFFER
              CBUFFER = "coordinates"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) NNodes
              WRITE (IO_UNIT%LEGO) (NODE(NPTUSED(n))%ID,                   n = 1,NNodes)
              WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Px,KIND(0E0)), n = 1,NNodes)
              WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Py,KIND(0E0)), n = 1,NNodes)
              WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Pz,KIND(0E0)), n = 1,NNodes)

              Nend = 0

              IF (MHX .GT. 0) THEN
                Nbgn = Nend + 1
                Nend = Nend + MHX
                CBUFFER = "hexa8"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) MHX
                WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
                WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
              ENDIF

              IF (MPX .GT. 0) THEN
                Nbgn = Nend + 1
                Nend = Nend + MPX
                CBUFFER = "penta6"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) MPX
                WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
                WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
              ENDIF

              IF (MPY .GT. 0) THEN
                Nbgn = Nend + 1
                Nend = Nend + MPY   
                CBUFFER = "pyramid5"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) MPY  
                WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
                WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
              ENDIF

              IF (MTX .GT. 0) THEN
                Nbgn = Nend + 1
                Nend = Nend + MTX   
                CBUFFER = "tetra4"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) MTX  
                WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
                WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
              ENDIF

              IF (MM3 .GT. 0) THEN
                Nbgn = Nend + 1
                Nend = Nend + MM3   
                CBUFFER = "tria3"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) MM3  
                WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
                WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
              ENDIF

              IF (MP3 .GT. 0) THEN
                Nbgn = Nend + 1
                Nend = Nend + MP3   
                CBUFFER = "tria3"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) MP3  
                WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
                WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
              ENDIF

              IF (MM4 .GT. 0) THEN
                Nbgn = Nend + 1
                Nend = Nend + MM4   
                CBUFFER = "quad4"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) MM4  
                WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
                WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
              ENDIF

              IF (MP4 .GT. 0) THEN
                Nbgn = Nend + 1
                Nend = Nend + MP4 
                CBUFFER = "quad4"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) MP4  
                WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
                WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
              ENDIF

              IF (MTR .GT. 0) THEN
                Nbgn = Nend + 1
                Nend = Nend + MTR   
                CBUFFER = "bar2"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) MTR  
                WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
                WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
              ENDIF

            ENDIF
          ENDDO
        ENDIF
!!
!! =========SECOND============================================================
!! Write the individual node-sets segment-sets, and wedge-sets as separate 
!! parts. The idea is that EnSight and ParaView can load these parts one at 
!! a time "on top of" the above mesh file to confirm that the individual sets 
!! have been correctly specified. (The underlying mesh should first be colored 
!! "gray" so that the colors assigned by ParaView to each set standout and are 
!! more easily examined for correctness.)
!!
        DO M = 1,NUMNS
!!
!! Clear marker/index-sequence and translation arrays.
!!
          NELUSED = 0
          NPTUSED = 0
          NPTNOWI = 0
!!
!! Mark nodes used by this node set.
!!
          N = 0
          DO WHILE (NEXT_NP_ID(M,N))
            NPTUSED(N) = 1
          ENDDO
!!
!! Covert NPTUSED (and NELUSED) into a sequential index map.
!!
          K = 0
          DO N = 1,NUMNP
            IF (NPTUSED(N) .EQ. 1) THEN
              K = K + 1
              NPTNOWI(N) = K
              NPTUSED(K) = N
              NELUSED(K) = N
            ENDIF
          ENDDO
!!
!! For later use, record the number of nodes and elements (each node will be 
!! a "vertex element" for ParaView) that will be in the file for this NP_SET.
!!
          NNodes = K
          NElems = K

          IF (NNodes .GT. 0) THEN
            WRITE (MEGO,'(I8.8)') NODE_SET(M)%SetID
!!
!! Node Set: Nodal point ID's and coordinates.
!!
            CBUFFER = "part"
            WRITE (IO_UNIT%LEGO) CBUFFER
            MPart = Mpart + 1
            WRITE (IO_UNIT%LEGO) MPart
            CBUFFER = "Node Set ID: "//MEGO
            WRITE (IO_UNIT%LEGO) CBUFFER
            CBUFFER = "coordinates"
            WRITE (IO_UNIT%LEGO) CBUFFER
            WRITE (IO_UNIT%LEGO) NNodes
            WRITE (IO_UNIT%LEGO) (NODE(NPTUSED(n))%ID,                   n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Px,KIND(0E0)), n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Py,KIND(0E0)), n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Pz,KIND(0E0)), n = 1,NNodes)
!!
!! Node Set: Point-element inventory.
!!
            CBUFFER = "point"
            WRITE (IO_UNIT%LEGO) CBUFFER
            WRITE (IO_UNIT%LEGO) NElems  
            WRITE (IO_UNIT%LEGO) (NODE(NPTUSED(n))%ID, n = 1,NElems)
            WRITE (IO_UNIT%LEGO) (n,                   n = 1,NElems)
  
          ENDIF
        ENDDO
!!
!! Now, write segment set files.
!!
        DO M = 1,NUMSS
!!
!! The function SEGMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! segments-used array NELUSED and esg/esg output-index translation array NPTNOWI for 
!! this segment set. If segment set M is not empty, the function returns 1, else 0. 
!!
          IF (SEGMENTS_AND_NODES_USED(M) .GT. 0) THEN
            WRITE (MEGO,'(I8.8)') SEGMENT_SET(M)%SetID
!!
!! Segment Set: Nodal point ID's and coordinates.
!!
            CBUFFER = "part"
            WRITE (IO_UNIT%LEGO) CBUFFER
            MPart = Mpart + 1
            WRITE (IO_UNIT%LEGO) MPart
            CBUFFER = "Segment Set ID: "//MEGO
            WRITE (IO_UNIT%LEGO) CBUFFER
            CBUFFER = "coordinates"
            WRITE (IO_UNIT%LEGO) CBUFFER
            WRITE (IO_UNIT%LEGO) NNodes
            WRITE (IO_UNIT%LEGO) (NODE(NPTUSED(n))%ID,                   n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Px,KIND(0E0)), n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Py,KIND(0E0)), n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Pz,KIND(0E0)), n = 1,NNodes)
!!
!! Segment Set: Quadrilateral and Triangle facet (element) inventory.
!!
            CBUFFER = "quad4"
            WRITE (IO_UNIT%LEGO) CBUFFER
            WRITE (IO_UNIT%LEGO) NElems  
            WRITE (IO_UNIT%LEGO) (SEGMENT(NELUSED(n))%PAR%SegID,         n = 1,NElems)
            WRITE (IO_UNIT%LEGO) (NTUPLE(1:SEGMENT_N_TUPLE(NELUSED(n))), n = 1,NElems)
  
          ENDIF
        ENDDO
!!
!! Now, write wedge-set files. (Wedge-sets are generated internally based
!! on each specified solid-to-solid tied interface (Keyword: MPCON4).
!!
        WedgeID = 0
        DO M = 1,NUMC4
!!
!! The function WEDGES_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED,
!! wedges-used array NELUSED and esg output-index translation array NPTNOWI for
!! this wedge set. If wedge set M is not empty, the function returns 1, else 0.
!!
          IF (WEDGES_AND_NODES_USED(M) .GT. 0) THEN
            WRITE (MEGO,'(I8.8)') SOLID_SOLID_INTERFACE(M)%MPCID
!!
!! Wedge Set: Nodal point ID's and coordinates.
!!
            CBUFFER = "part"
            WRITE (IO_UNIT%LEGO) CBUFFER
            MPart = Mpart + 1
            WRITE (IO_UNIT%LEGO) MPart
            CBUFFER = "Wedge Set ID: "//MEGO
            WRITE (IO_UNIT%LEGO) CBUFFER
            CBUFFER = "coordinates"
            WRITE (IO_UNIT%LEGO) CBUFFER
            WRITE (IO_UNIT%LEGO) NNodes
            WRITE (IO_UNIT%LEGO) (n,                                        n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(WEDGE_NP_COORDINATE(n,1),KIND(0E0)), n = 1,NNodes)  !  x-coordinate
            WRITE (IO_UNIT%LEGO) (REAL(WEDGE_NP_COORDINATE(n,2),KIND(0E0)), n = 1,NNodes)  !  y-coordinate
            WRITE (IO_UNIT%LEGO) (REAL(WEDGE_NP_COORDINATE(n,3),KIND(0E0)), n = 1,NNodes)  !  z-coordinate
!!
!! Wedge Set: Interface wedge set inventory for this tied-contact (MPCON4).
!!
            CBUFFER = "penta6"
            WRITE (IO_UNIT%LEGO) CBUFFER
            WRITE (IO_UNIT%LEGO) NElems  
            WRITE (IO_UNIT%LEGO) (WedgeID + n,                n = 1,NElems)
            WRITE (IO_UNIT%LEGO) (NTUPLE(1:WEDGE_N_TUPLE(n)), n = 1,NElems)

            WedgeID = WedgeID + NElems

          ENDIF
        ENDDO

        CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')
!!
!! =========THIRD=============================================================
!! Write nodal point and cell (element) data.
!! 
!! Open an ESG nodal point velocity initial condition file.
!!
        IOERROR = .TRUE.
        OPEN
     &    (
     &    UNIT    =  IO_UNIT%LEGO,
     &    FILE    = 'fmaego.mesh.vics',
     &    STATUS  = 'UNKNOWN',
#ifdef _G95_
     &    FORM    = 'UNFORMATTED', ACCESS='STREAM',  !  G95
#endif
#ifdef _CVF_NT_
     &    FORM    = 'BINARY', CONVERT='BIG_ENDIAN',  !  CVF
#endif
#ifdef LANGUAGE_FORTRAN90
     &    FORM    = 'BINARY',                        !  Pathf95
#endif
     &    ERR     =  200
     &    )
        IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 200    IF (IOERROR) THEN
          CALL USER_MESSAGE
     &      (
     &      MSGL//'FATAL'//
     &      MSGL//'WRITE_ESG_DATA_FILE.002.00'//
     &      MSGL//'Unable To Execute OPEN On: '//'fmaego.mesh.vics'
     &      )
        ELSE
!!
!! Initialize sequential Material part counter.
!!
          MPart = 0
!!
!! Initialize EnSight Gold static variable results file
!!
          CBUFFER = "Velocity Initial Conditions"
          WRITE (IO_UNIT%LEGO) CBUFFER
!!
!! Loop on Material parts.
!!
          DO M = 1,NUMMT
!!
!! The function ELEMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! elements-used array NELUSED and vtk output-index translation array NPTNOWI for 
!! this material. If material M is used, the function returns 1, otherwise 0. 
!!
            IF (ELEMENTS_AND_NODES_USED(M) .GT. 0) THEN

              NNodes = MATERIAL(M)%NNodes
              NElems = MATERIAL(M)%NElems

              CBUFFER = "part"
              WRITE (IO_UNIT%LEGO) CBUFFER
              MPart = MPart + 1
              WRITE (IO_UNIT%LEGO) MPart
              CBUFFER = "coordinates"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Vx,KIND(0E0)), n = 1,NNodes)
              WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Vy,KIND(0E0)), n = 1,NNodes)
              WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Vz,KIND(0E0)), n = 1,NNodes)

            ENDIF
          ENDDO
!!
!! Assign initial velocity values to node sets.
!!
          DO M = 1,NUMNS
!!
!! Clear marker/index-sequence and translation arrays.
!!
            NELUSED = 0
            NPTUSED = 0
            NPTNOWI = 0
!!
!! Mark nodes used by this node set.
!!
            N = 0
            DO WHILE (NEXT_NP_ID(M,N))
              NPTUSED(N) = 1
            ENDDO
!!
!! Covert NPTUSED (and NELUSED) into a sequential index map.
!!
            K = 0
            DO N = 1,NUMNP
              IF (NPTUSED(N) .EQ. 1) THEN
                K = K + 1
                NPTNOWI(N) = K
                NPTUSED(K) = N
                NELUSED(K) = N
              ENDIF
            ENDDO
!!
!! For later use, record the number of nodes and elements (each node will be 
!! a "vertex element" for ParaView) that will be in the file for this NP_SET.
!!
            NNodes = K
            NElems = K

            IF (NNodes .GT. 0) THEN
!!
!! Node Set: Velocity initial conditions.
!!
              CBUFFER = "part"
              WRITE (IO_UNIT%LEGO) CBUFFER
              MPart = MPart + 1
              WRITE (IO_UNIT%LEGO) MPart
              CBUFFER = "coordinates"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Vx,KIND(0E0)), n = 1,NNodes)
              WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Vy,KIND(0E0)), n = 1,NNodes)
              WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Vz,KIND(0E0)), n = 1,NNodes)

            ENDIF
          ENDDO
!!
!! Assign initial velocity values to segment sets.
!!
          DO M = 1,NUMSS
!!
!! The function SEGMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! segments-used array NELUSED and esg/esg output-index translation array NPTNOWI for 
!! this segment set. If segment set M is not empty, the function returns 1, else 0. 
!!
            IF (SEGMENTS_AND_NODES_USED(M) .GT. 0) THEN
!!
!! Segment Set: Velocity initial conditions.
!!
              CBUFFER = "part"
              WRITE (IO_UNIT%LEGO) CBUFFER
              MPart = Mpart + 1
              WRITE (IO_UNIT%LEGO) MPart
              CBUFFER = "coordinates"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Vx,KIND(0E0)), n = 1,NNodes)
              WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Vy,KIND(0E0)), n = 1,NNodes)
              WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Vz,KIND(0E0)), n = 1,NNodes)

            ENDIF
          ENDDO
!!
!! Assign initial velocity values to wedge sets.
!!
          DO M = 1,NUMC4
!!
!! The function WEDGES_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED,
!! wedges-used array NELUSED and esg output-index translation array NPTNOWI for
!! this wedge set. If wedge set M is not empty, the function returns 1, else 0.
!!
            IF (WEDGES_AND_NODES_USED(M) .GT. 0) THEN
!!
!! Wedge Set: Velocity initial conditions.
!!
              CBUFFER = "part"
              WRITE (IO_UNIT%LEGO) CBUFFER
              MPart = Mpart + 1
              WRITE (IO_UNIT%LEGO) MPart
              CBUFFER = "coordinates"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(WEDGE_NP_VELOCITY(n,1),KIND(0E0)), n = 1,NNodes)  !  x-coordinate
              WRITE (IO_UNIT%LEGO) (REAL(WEDGE_NP_VELOCITY(n,2),KIND(0E0)), n = 1,NNodes)  !  y-coordinate
              WRITE (IO_UNIT%LEGO) (REAL(WEDGE_NP_VELOCITY(n,3),KIND(0E0)), n = 1,NNodes)  !  z-coordinate

            ENDIF
          ENDDO

          CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')
        ENDIF
!! 
!! Open an ESG element material number file.
!!
        IOERROR = .TRUE.
        OPEN
     &    (
     &    UNIT    =  IO_UNIT%LEGO,
     &    FILE    = 'fmaego.mesh.mats',
     &    STATUS  = 'UNKNOWN',
#ifdef _G95_
     &    FORM    = 'UNFORMATTED', ACCESS='STREAM',  !  G95
#endif
#ifdef _CVF_NT_
     &    FORM    = 'BINARY', CONVERT='BIG_ENDIAN',  !  CVF
#endif
#ifdef LANGUAGE_FORTRAN90
     &    FORM    = 'BINARY',                        !  Pathf95
#endif
     &    ERR     =  300
     &    )
        IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 300    IF (IOERROR) THEN
          CALL USER_MESSAGE
     &      (
     &      MSGL//'FATAL'//
     &      MSGL//'WRITE_ESG_DATA_FILE.003.00'//
     &      MSGL//'Unable To Execute OPEN On: '//'fmaego.mesh.mats'
     &      )
        ELSE
!!
!! Initialize sequential part counter.
!!
          MPart = 0
!!
!! Initialize EnSight Gold static variable results file.
!!
          CBUFFER = "Cell Material Color Index"
          WRITE (IO_UNIT%LEGO) CBUFFER
!!
!! Loop on Material parts.
!!
          DO M = 1,NUMMT
!!
!! The function ELEMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! elements-used array NELUSED and vtk output-index translation array NPTNOWI for 
!! this material. If material M is used, the function returns 1, otherwise 0. 
!!
            IF (ELEMENTS_AND_NODES_USED(M) .GT. 0) THEN

              NNodes = MATERIAL(M)%NNodes
              NElems = MATERIAL(M)%NElems

              CBUFFER = "part"
              WRITE (IO_UNIT%LEGO) CBUFFER
              MPart = MPart + 1
              WRITE (IO_UNIT%LEGO) MPart

              IF (MHX .GT. 0) THEN
                CBUFFER = "hexa8"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) (REAL(M,KIND(0E0)), n = 1,MHX)
              ENDIF

              IF (MPX .GT. 0) THEN
                CBUFFER = "penta6"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) (REAL(M,KIND(0E0)), n = 1,MPX)
              ENDIF

              IF (MPY .GT. 0) THEN
                CBUFFER = "pyramid5"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) (REAL(M,KIND(0E0)), n = 1,MPY)
              ENDIF

              IF (MTX .GT. 0) THEN
                CBUFFER = "tetra4"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) (REAL(M,KIND(0E0)), n = 1,MTX)
              ENDIF

              IF (MM3 .GT. 0) THEN
                CBUFFER = "tria3"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) (REAL(M,KIND(0E0)), n = 1,MM3)
              ENDIF

              IF (MP3 .GT. 0) THEN
                CBUFFER = "tria3"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) (REAL(M,KIND(0E0)), n = 1,MP3)
              ENDIF

              IF (MM4 .GT. 0) THEN
                CBUFFER = "quad4"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) (REAL(M,KIND(0E0)), n = 1,MM4)
              ENDIF

              IF (MP4 .GT. 0) THEN
                CBUFFER = "quad4"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) (REAL(M,KIND(0E0)), n = 1,MP4)
              ENDIF

              IF (MTR .GT. 0) THEN
                CBUFFER = "bar2"
                WRITE (IO_UNIT%LEGO) CBUFFER
                WRITE (IO_UNIT%LEGO) (REAL(M,KIND(0E0)), n = 1,MTR)
              ENDIF

            ENDIF
          ENDDO
!!
!! Assign "MatID" values to node sets.
!!
          DO M = 1,NUMNS
!!
!! Clear marker/index-sequence and translation arrays.
!!
            NELUSED = 0
            NPTUSED = 0
            NPTNOWI = 0
!!
!! Mark nodes used by this node set.
!!
            N = 0
            DO WHILE (NEXT_NP_ID(M,N))
              NPTUSED(N) = 1
            ENDDO
!!
!! Covert NPTUSED (and NELUSED) into a sequential index map.
!!
            K = 0
            DO N = 1,NUMNP
              IF (NPTUSED(N) .EQ. 1) THEN
                K = K + 1
                NPTNOWI(N) = K
                NPTUSED(K) = N
                NELUSED(K) = N
              ENDIF
            ENDDO
!!
!! For later use, record the number of nodes and elements (each node will be 
!! a "vertex element" for ParaView) that will be in the file for this NP_SET.
!!
            NNodes = K
            NElems = K

            IF (NNodes .GT. 0) THEN

              CBUFFER = "part"
              WRITE (IO_UNIT%LEGO) CBUFFER
              MPart = Mpart + 1
              WRITE (IO_UNIT%LEGO) MPart
!!
!! Node Set: Assign internal node set ID as "MatID."
!!
              CBUFFER = "point"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(NUMMT+M,KIND(0E0)), n = 1,NElems)

            ENDIF
          ENDDO
!!
!! Assign "MatID" values to segment sets.
!!
          DO M = 1,NUMSS
!!
!! The function SEGMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! segments-used array NELUSED and esg/esg output-index translation array NPTNOWI for 
!! this segment set. If segment set M is not empty, the function returns 1, else 0. 
!!
            IF (SEGMENTS_AND_NODES_USED(M) .GT. 0) THEN

              CBUFFER = "part"
              WRITE (IO_UNIT%LEGO) CBUFFER
              MPart = Mpart + 1
              WRITE (IO_UNIT%LEGO) MPart
!!
!! Segment Set: Assign internal segment set ID as "MatID."
!!
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(NUMMT+NUMNS+M,KIND(0E0)), n = 1,NElems)

            ENDIF
          ENDDO
!!
!! Assign "MatID" values to wedge sets.
!!
          DO M = 1,NUMC4
!!
!! The function WEDGES_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED,
!! wedges-used array NELUSED and esg output-index translation array NPTNOWI for
!! this wedge set. If wedge set M is not empty, the function returns 1, else 0.
!!
            IF (WEDGES_AND_NODES_USED(M) .GT. 0) THEN

              CBUFFER = "part"
              WRITE (IO_UNIT%LEGO) CBUFFER
              MPart = Mpart + 1
              WRITE (IO_UNIT%LEGO) MPart
!!
!! Wedge Set: Assign internal wedge set ID as "MatID."
!!
              CBUFFER = "penta6"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(NUMMT+NUMNS+NUMSS+M,KIND(0E0)), n = 1,NElems)

            ENDIF
          ENDDO

          CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')
        ENDIF
!!
!! =========FOURTH============================================================
!! Write a *.case file to allow EnSight and ParaView to identify the mesh, 
!! node-set, side-set and wedge-set files as distinct "Parts."
!! 
!! Open an ESG case-file.
!!
        IOERROR = .TRUE.
        OPEN
     &    (
     &    UNIT    =  IO_UNIT%LEGO,
     &    FILE    = 'fmaego.mesh.case',
     &    STATUS  = 'UNKNOWN',
     &    FORM    = 'FORMATTED',
     &    ERR     =  400
     &    )
        IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 400    IF (IOERROR) THEN
          CALL USER_MESSAGE
     &      (
     &      MSGL//'FATAL'//
     &      MSGL//'WRITE_ESG_DATA_FILE.004.00'//
     &      MSGL//'Unable To Execute OPEN On: '//'fmaego.mesh.case'
     &      )
        ELSE

          WRITE (IO_UNIT%LEGO,'(A)') "FORMAT"
          WRITE (IO_UNIT%LEGO,'(A)') "type: ensight gold"

          WRITE (IO_UNIT%LEGO,'(A)') "GEOMETRY"
          WRITE (IO_UNIT%LEGO,'(A)') "model: fmaego.mesh.geom"

          WRITE (IO_UNIT%LEGO,'(A)') "VARIABLE"
          WRITE (IO_UNIT%LEGO,'(A)') "vector per node: Velocity-IC fmaego.mesh.vics"
          WRITE (IO_UNIT%LEGO,'(A)') "scalar per element: Material fmaego.mesh.mats"

          CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')
        ENDIF

        FIRST = .FALSE.
      ENDIF
!!
!!
!!############################################################################
!! PART 2. RESULTS FOR THIS TIME STEP.
!! Initialize and append nodal point and cell (element) data.
!!
!! Increment EnSight Gold number-of-time-steps counter and store current time.
!!
      TStep = TStep + 1
      IF (TStep .LE. MStep) THEN
        WRITE (NEGO,'(I5)') TStep
        TIME_VALUES(TStep) = TIME%Total
      ELSE
        WRITE (MSG1,'(I8)') TStep
        WRITE (MSG2,'(I8)') MStep
        CALL USER_MESSAGE
     &    (
     &    MSGL//'FATAL'//
     &    MSGL//'WRITE_ESG_DATA_FILE.010.02'//
     &    MSGL//'TStep, The Current Write Counter, Has Exceeded MStep.'//
     &    MSGL//'TStep:'//MSG1//
     &    MSGL//'MStep:'//MSG2
     &    )
      ENDIF
!!
!! =========FIRST=============================================================
!! Open geometry file fmaego.data.geom (EnSight Gold output) for nodal point 
!! and cell results to reference.
!!
      IOERROR = .TRUE.
      OPEN
     &  (
     &  UNIT     =  IO_UNIT%LEGO,
     &  FILE     = 'fmaego.data.geom',
     &  STATUS   = 'UNKNOWN',
#ifdef _G95_
     &  FORM     = 'UNFORMATTED', ACCESS='STREAM',  !  G95
#endif
#ifdef _CVF_NT_
     &  FORM     = 'BINARY', CONVERT='BIG_ENDIAN',  !  CVF
#endif
#ifdef LANGUAGE_FORTRAN90
     &  FORM     = 'BINARY',                        !  Pathf95
#endif
     &  POSITION = 'APPEND',
     &  ERR      =  500
     &  )
      IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 500  IF (IOERROR) THEN
        CALL USER_MESSAGE
     &    (
     &    MSGL//'FATAL'//
     &    MSGL//'WRITE_ESG_DATA_FILE.005.00'//
     &    MSGL//'Unable To Execute OPEN On: '//'fmaego.data.geom'
     &    )
      ELSE
!!
!! Initialize sequential Material part counter.
!!
        MPart = 0
!!
!! Insert EnSight Gold C-binary data file que (see FORM='UNFORMATTED' above).
!!
        IF (INSERT_C_BINARY_QUE) THEN
          CBUFFER = "C Binary"                       ! EnSight reader format que.
          WRITE (IO_UNIT%LEGO) CBUFFER
          INSERT_C_BINARY_QUE = .FALSE.
        ENDIF
!!
!! Start this-time-step block.
!!
        CBUFFER = "BEGIN TIME STEP"                  ! "Time serialized" que.
        WRITE (IO_UNIT%LEGO) CBUFFER
        CBUFFER = TRIM(JOB_ID_RECORD%CURRENT%TITLE)  ! User's job title record.
        WRITE (IO_UNIT%LEGO) CBUFFER
        CBUFFER = " "                                ! (unused subtitle)
        WRITE (IO_UNIT%LEGO) CBUFFER
        CBUFFER = "node id given"                    ! Expect to read nodal ID's
        WRITE (IO_UNIT%LEGO) CBUFFER
        CBUFFER = "element id given"                 ! Expect to read element ID's
        WRITE (IO_UNIT%LEGO) CBUFFER
!!
!! Turn individual material domains into esg/vtk "parts."
!!
        DO M = 1,NUMMT
!!
!! The function ELEMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! elements-used array NELUSED and vtk output-index translation array NPTNOWI for 
!! this material. If material M is used, the function returns 1, otherwise 0. 
!!
          IF (ELEMENTS_AND_NODES_USED(M) .GT. 0) THEN

            WRITE (MEGO,'(I8.8)') MATERIAL(M)%MatID

            NNodes = MATERIAL(M)%NNodes
            NElems = MATERIAL(M)%NElems

            CBUFFER = "part"
            WRITE (IO_UNIT%LEGO) CBUFFER
            MPart = Mpart + 1
            WRITE (IO_UNIT%LEGO) MPart
            CBUFFER = "Part with MatID: "//MEGO
            WRITE (IO_UNIT%LEGO) CBUFFER
            CBUFFER = "coordinates"
            WRITE (IO_UNIT%LEGO) CBUFFER
            WRITE (IO_UNIT%LEGO) NNodes
            WRITE (IO_UNIT%LEGO) (NODE(NPTUSED(n))%ID,                                         n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Px+MOTION(NPTUSED(n))%Ux,KIND(0E0)), n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Py+MOTION(NPTUSED(n))%Uy,KIND(0E0)), n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Pz+MOTION(NPTUSED(n))%Uz,KIND(0E0)), n = 1,NNodes)

            Nend = 0

            IF (MHX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MHX
              CBUFFER = "hexa8"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) MHX
              WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
            ENDIF

            IF (MPX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MPX
              CBUFFER = "penta6"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) MPX
              WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
            ENDIF

            IF (MPY .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MPY   
              CBUFFER = "pyramid5"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) MPY  
              WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
            ENDIF

            IF (MTX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MTX   
              CBUFFER = "tetra4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) MTX  
              WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
            ENDIF

            IF (MM3 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MM3   
              CBUFFER = "tria3"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) MM3  
              WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
            ENDIF

            IF (MP3 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MP3   
              CBUFFER = "tria3"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) MP3  
              WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
            ENDIF

            IF (MM4 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MM4   
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) MM4  
              WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
            ENDIF

            IF (MP4 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MP4 
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) MP4  
              WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
            ENDIF

            IF (MTR .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MTR   
              CBUFFER = "bar2"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) MTR  
              WRITE (IO_UNIT%LEGO) (EleID_DATA(NELUSED(n)),                n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (NTUPLE(1:ELEMENT_N_TUPLE(NELUSED(n))), n = Nbgn,Nend)
            ENDIF

          ENDIF
        ENDDO
!!
!! Close out this-time-step block.
!!
        CBUFFER = "END TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')
      ENDIF
!!
!! =========SECOND============================================================
!! Open and append ESG nodal point results to data files.
!!
!! Displacements...
!!
      IOERROR = .TRUE.
      OPEN
     &  (
     &  UNIT     =  IO_UNIT%LEGO,
     &  FILE     = 'fmaego.data.ndis',
     &  STATUS   = 'UNKNOWN',
#ifdef _G95_
     &  FORM     = 'UNFORMATTED', ACCESS='STREAM',  !  G95
#endif
#ifdef _CVF_NT_
     &  FORM     = 'BINARY', CONVERT='BIG_ENDIAN',  !  CVF
#endif
#ifdef LANGUAGE_FORTRAN90
     &  FORM     = 'BINARY',                        !  Pathf95
#endif
     &  POSITION = 'APPEND',
     &  ERR      =  601
     &  )
      IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 601  IF (IOERROR) THEN
        CALL USER_MESSAGE
     &    (
     &    MSGL//'FATAL'//
     &    MSGL//'WRITE_ESG_DATA_FILE.006.01'//
     &    MSGL//'Unable To Execute OPEN On: '//'fmaego.data.ndis'
     &    )
      ELSE
!!
!! Initialize sequential Material part counter.
!!
        MPart = 0
!!
!! Start this-time-step block.
!!
        CBUFFER = "BEGIN TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CBUFFER = "Nodal Point Displacement Results"
        WRITE (IO_UNIT%LEGO) CBUFFER
!!
!! Loop on Material parts.
!!
        DO M = 1,NUMMT
!!
!! The function ELEMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! elements-used array NELUSED and vtk output-index translation array NPTNOWI for 
!! this material. If material M is used, the function returns 1, otherwise 0. 
!!
          IF (ELEMENTS_AND_NODES_USED(M) .GT. 0) THEN

            NNodes = MATERIAL(M)%NNodes
            NElems = MATERIAL(M)%NElems

            CBUFFER = "part"
            WRITE (IO_UNIT%LEGO) CBUFFER
            MPart = MPart + 1
            WRITE (IO_UNIT%LEGO) MPart
            CBUFFER = "coordinates"
            WRITE (IO_UNIT%LEGO) CBUFFER
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Ux,KIND(0E0)), n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Uy,KIND(0E0)), n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Uz,KIND(0E0)), n = 1,NNodes)

          ENDIF
        ENDDO
!!
!! Close out this-time-step block.
!!
        CBUFFER = "END TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')

      ENDIF
!!
!! Velocities...
!!
      IOERROR = .TRUE.
      OPEN
     &  (
     &  UNIT     =  IO_UNIT%LEGO,
     &  FILE     = 'fmaego.data.nvel',
     &  STATUS   = 'UNKNOWN',
#ifdef _G95_
     &  FORM     = 'UNFORMATTED', ACCESS='STREAM',  !  G95
#endif
#ifdef _CVF_NT_
     &  FORM     = 'BINARY', CONVERT='BIG_ENDIAN',  !  CVF
#endif
#ifdef LANGUAGE_FORTRAN90
     &  FORM     = 'BINARY',                        !  Pathf95
#endif
     &  POSITION = 'APPEND',
     &  ERR      =  602
     &  )
      IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 602  IF (IOERROR) THEN
        CALL USER_MESSAGE
     &    (
     &    MSGL//'FATAL'//
     &    MSGL//'WRITE_ESG_DATA_FILE.006.02'//
     &    MSGL//'Unable To Execute OPEN On: '//'fmaego.data.nvel'
     &    )
      ELSE
!!
!! Initialize sequential Material part counter.
!!
        MPart = 0
!!
!! Start this-time-step block.
!!
        CBUFFER = "BEGIN TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CBUFFER = "Nodal Point Displacement Results"
        WRITE (IO_UNIT%LEGO) CBUFFER
!!
!! Loop on Material parts.
!!
        DO M = 1,NUMMT
!!
!! The function ELEMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! elements-used array NELUSED and vtk output-index translation array NPTNOWI for 
!! this material. If material M is used, the function returns 1, otherwise 0. 
!!
          IF (ELEMENTS_AND_NODES_USED(M) .GT. 0) THEN

            NNodes = MATERIAL(M)%NNodes
            NElems = MATERIAL(M)%NElems

            CBUFFER = "part"
            WRITE (IO_UNIT%LEGO) CBUFFER
            MPart = MPart + 1
            WRITE (IO_UNIT%LEGO) MPart
            CBUFFER = "coordinates"
            WRITE (IO_UNIT%LEGO) CBUFFER
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Vx,KIND(0E0)), n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Vy,KIND(0E0)), n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Vz,KIND(0E0)), n = 1,NNodes)

          ENDIF
        ENDDO
!!
!! Close out this-time-step block.
!!
        CBUFFER = "END TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')

      ENDIF
!!
!! Accelerations...
!!
      IOERROR = .TRUE.
      OPEN
     &  (
     &  UNIT     =  IO_UNIT%LEGO,
     &  FILE     = 'fmaego.data.nacc',
     &  STATUS   = 'UNKNOWN',
#ifdef _G95_
     &  FORM     = 'UNFORMATTED', ACCESS='STREAM',  !  G95
#endif
#ifdef _CVF_NT_
     &  FORM     = 'BINARY', CONVERT='BIG_ENDIAN',  !  CVF
#endif
#ifdef LANGUAGE_FORTRAN90
     &  FORM     = 'BINARY',                        !  Pathf95
#endif
     &  POSITION = 'APPEND',
     &  ERR      =  603
     &  )
      IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 603  IF (IOERROR) THEN
        CALL USER_MESSAGE
     &    (
     &    MSGL//'FATAL'//
     &    MSGL//'WRITE_ESG_DATA_FILE.006.03'//
     &    MSGL//'Unable To Execute OPEN On: '//'fmaego.data.nacc'
     &    )
      ELSE
!!
!! Initialize sequential Material part counter.
!!
        MPart = 0
!!
!! Start this-time-step block.
!!
        CBUFFER = "BEGIN TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CBUFFER = "Nodal Point Displacement Results"
        WRITE (IO_UNIT%LEGO) CBUFFER
!!
!! Loop on Material parts.
!!
        DO M = 1,NUMMT
!!
!! The function ELEMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! elements-used array NELUSED and vtk output-index translation array NPTNOWI for 
!! this material. If material M is used, the function returns 1, otherwise 0. 
!!
          IF (ELEMENTS_AND_NODES_USED(M) .GT. 0) THEN

            NNodes = MATERIAL(M)%NNodes
            NElems = MATERIAL(M)%NElems

            CBUFFER = "part"
            WRITE (IO_UNIT%LEGO) CBUFFER
            MPart = MPart + 1
            WRITE (IO_UNIT%LEGO) MPart
            CBUFFER = "coordinates"
            WRITE (IO_UNIT%LEGO) CBUFFER
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Ax,KIND(0E0)), n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Ay,KIND(0E0)), n = 1,NNodes)
            WRITE (IO_UNIT%LEGO) (REAL(MOTION(NPTUSED(n))%Az,KIND(0E0)), n = 1,NNodes)

          ENDIF
        ENDDO
!!
!! Close out this-time-step block.
!!
        CBUFFER = "END TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')

      ENDIF
!!
!! =========THIRD=============================================================
!! Open and append ESG element (cell) results to data files.
!!
!! Internal Material Number...
!!
      IOERROR = .TRUE.
      OPEN
     &  (
     &  UNIT     =  IO_UNIT%LEGO,
     &  FILE     = 'fmaego.data.emat',
     &  STATUS   = 'UNKNOWN',
#ifdef _G95_
     &  FORM     = 'UNFORMATTED', ACCESS='STREAM',  !  G95
#endif
#ifdef _CVF_NT_
     &  FORM     = 'BINARY', CONVERT='BIG_ENDIAN',  !  CVF
#endif
#ifdef LANGUAGE_FORTRAN90
     &  FORM     = 'BINARY',                        !  Pathf95
#endif
     &  POSITION = 'APPEND',
     &  ERR      =  701
     &  )
      IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 701  IF (IOERROR) THEN
        CALL USER_MESSAGE
     &    (
     &    MSGL//'FATAL'//
     &    MSGL//'WRITE_ESG_DATA_FILE.007.01'//
     &    MSGL//'Unable To Execute OPEN On: '//'fmaego.data.emat'
     &    )
      ELSE
!!
!! Initialize sequential Material part counter.
!!
        MPart = 0
!!
!! Start this-time-step block.
!!
        CBUFFER = "BEGIN TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CBUFFER = "Element Internal Material Index"
        WRITE (IO_UNIT%LEGO) CBUFFER
!!
!! Loop on Material parts.
!!
        DO M = 1,NUMMT
!!
!! The function ELEMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! elements-used array NELUSED and vtk output-index translation array NPTNOWI for 
!! this material. If material M is used, the function returns 1, otherwise 0. 
!!
          IF (ELEMENTS_AND_NODES_USED(M) .GT. 0) THEN

            NNodes = MATERIAL(M)%NNodes
            NElems = MATERIAL(M)%NElems

            CBUFFER = "part"
            WRITE (IO_UNIT%LEGO) CBUFFER
            MPart = MPart + 1
            WRITE (IO_UNIT%LEGO) MPart

            Nend = 0

            IF (MHX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MHX
              CBUFFER = "hexa8"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(MatID_DATA(NELUSED(n)),KIND(0E0)), n = Nbgn,Nend)
            ENDIF

            IF (MPX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MPX
              CBUFFER = "penta6"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(MatID_DATA(NELUSED(n)),KIND(0E0)), n = Nbgn,Nend)
            ENDIF

            IF (MPY .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MPY   
              CBUFFER = "pyramid5"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(MatID_DATA(NELUSED(n)),KIND(0E0)), n = Nbgn,Nend)
            ENDIF

            IF (MTX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MTX   
              CBUFFER = "tetra4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(MatID_DATA(NELUSED(n)),KIND(0E0)), n = Nbgn,Nend)
            ENDIF

            IF (MM3 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MM3   
              CBUFFER = "tria3"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(MatID_DATA(NELUSED(n)),KIND(0E0)), n = Nbgn,Nend)
            ENDIF

            IF (MP3 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MP3   
              CBUFFER = "tria3"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(MatID_DATA(NELUSED(n)),KIND(0E0)), n = Nbgn,Nend)
            ENDIF

            IF (MM4 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MM4   
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(MatID_DATA(NELUSED(n)),KIND(0E0)), n = Nbgn,Nend)
            ENDIF

            IF (MP4 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MP4 
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(MatID_DATA(NELUSED(n)),KIND(0E0)), n = Nbgn,Nend)
            ENDIF

            IF (MTR .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MTR   
              CBUFFER = "bar2"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (REAL(MatID_DATA(NELUSED(n)),KIND(0E0)), n = Nbgn,Nend)
            ENDIF

          ENDIF
        ENDDO
!!
!! Close out this-time-step block.
!!
        CBUFFER = "END TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')

      ENDIF
!!
!! Stress...
!!
      IOERROR = .TRUE.
      OPEN
     &  (
     &  UNIT     =  IO_UNIT%LEGO,
     &  FILE     = 'fmaego.data.estr',
     &  STATUS   = 'UNKNOWN',
#ifdef _G95_
     &  FORM     = 'UNFORMATTED', ACCESS='STREAM',  !  G95
#endif
#ifdef _CVF_NT_
     &  FORM     = 'BINARY', CONVERT='BIG_ENDIAN',  !  CVF
#endif
#ifdef LANGUAGE_FORTRAN90
     &  FORM     = 'BINARY',                        !  Pathf95
#endif
     &  POSITION = 'APPEND',
     &  ERR      =  702
     &  )
      IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 702  IF (IOERROR) THEN
        CALL USER_MESSAGE
     &    (
     &    MSGL//'FATAL'//
     &    MSGL//'WRITE_ESG_DATA_FILE.007.02'//
     &    MSGL//'Unable To Execute OPEN On: '//'fmaego.data.estr'
     &    )
      ELSE
!!
!! Initialize sequential Material part counter.
!!
        MPart = 0
!!
!! Start this-time-step block.
!!
        CBUFFER = "BEGIN TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CBUFFER = "Element Stress(xx,yy,zz,xy,xz,yz) Results"
        WRITE (IO_UNIT%LEGO) CBUFFER
!!
!! Loop on Material parts.
!!
        DO M = 1,NUMMT
!!
!! The function ELEMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! elements-used array NELUSED and vtk output-index translation array NPTNOWI for 
!! this material. If material M is used, the function returns 1, otherwise 0. 
!!
          IF (ELEMENTS_AND_NODES_USED(M) .GT. 0) THEN

            NNodes = MATERIAL(M)%NNodes
            NElems = MATERIAL(M)%NElems

            CBUFFER = "part"
            WRITE (IO_UNIT%LEGO) CBUFFER
            MPart = MPart + 1
            WRITE (IO_UNIT%LEGO) MPart

            Nend = 0

            IF (MHX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MHX
              CBUFFER = "hexa8"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),1), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),2), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),3), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),4), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),5), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),6), n = Nbgn,Nend)
            ENDIF

            IF (MPX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MPX
              CBUFFER = "penta6"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),1), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),2), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),3), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),4), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),5), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),6), n = Nbgn,Nend)
            ENDIF

            IF (MPY .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MPY   
              CBUFFER = "pyramid5"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),1), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),2), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),3), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),4), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),5), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),6), n = Nbgn,Nend)
            ENDIF

            IF (MTX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MTX   
              CBUFFER = "tetra4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),1), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),2), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),3), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),4), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),5), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),6), n = Nbgn,Nend)
            ENDIF

            IF (MM3 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MM3   
              CBUFFER = "tria3"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),1), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),2), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),3), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),4), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),5), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),6), n = Nbgn,Nend)
            ENDIF

            IF (MP3 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MP3   
              CBUFFER = "tria3"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),1), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),2), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),3), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),4), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),5), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),6), n = Nbgn,Nend)
            ENDIF

            IF (MM4 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MM4   
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),1), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),2), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),3), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),4), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),5), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),6), n = Nbgn,Nend)
            ENDIF

            IF (MP4 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MP4 
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),1), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),2), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),3), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),4), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),5), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),6), n = Nbgn,Nend)
            ENDIF

            IF (MTR .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MTR   
              CBUFFER = "bar2"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),1), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),2), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),3), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),4), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),5), n = Nbgn,Nend)
              WRITE (IO_UNIT%LEGO) (STRESS_DATA(NELUSED(n),6), n = Nbgn,Nend)
            ENDIF

          ENDIF
        ENDDO
!!
!! Close out this-time-step block.
!!
        CBUFFER = "END TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')

      ENDIF
!!
!! Bulk Strain...
!!
      IOERROR = .TRUE.
      OPEN
     &  (
     &  UNIT     =  IO_UNIT%LEGO,
     &  FILE     = 'fmaego.data.elnv',
     &  STATUS   = 'UNKNOWN',
#ifdef _G95_
     &  FORM     = 'UNFORMATTED', ACCESS='STREAM',  !  G95
#endif
#ifdef _CVF_NT_
     &  FORM     = 'BINARY', CONVERT='BIG_ENDIAN',  !  CVF
#endif
#ifdef LANGUAGE_FORTRAN90
     &  FORM     = 'BINARY',                        !  Pathf95
#endif
     &  POSITION = 'APPEND',
     &  ERR      =  703
     &  )
      IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 703  IF (IOERROR) THEN
        CALL USER_MESSAGE
     &    (
     &    MSGL//'FATAL'//
     &    MSGL//'WRITE_ESG_DATA_FILE.007.03'//
     &    MSGL//'Unable To Execute OPEN On: '//'fmaego.data.elnv'
     &    )
      ELSE
!!
!! Initialize sequential Material part counter.
!!
        MPart = 0
!!
!! Start this-time-step block.
!!
        CBUFFER = "BEGIN TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CBUFFER = "Element Bulk Strain Results"
        WRITE (IO_UNIT%LEGO) CBUFFER
!!
!! Loop on Material parts.
!!
        DO M = 1,NUMMT
!!
!! The function ELEMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! elements-used array NELUSED and vtk output-index translation array NPTNOWI for 
!! this material. If material M is used, the function returns 1, otherwise 0. 
!!
          IF (ELEMENTS_AND_NODES_USED(M) .GT. 0) THEN

            NNodes = MATERIAL(M)%NNodes
            NElems = MATERIAL(M)%NElems

            CBUFFER = "part"
            WRITE (IO_UNIT%LEGO) CBUFFER
            MPart = MPart + 1
            WRITE (IO_UNIT%LEGO) MPart

            Nend = 0

            IF (MHX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MHX
              CBUFFER = "hexa8"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (BULK_STRAIN(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MPX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MPX
              CBUFFER = "penta6"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (BULK_STRAIN(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MPY .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MPY   
              CBUFFER = "pyramid5"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (BULK_STRAIN(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MTX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MTX   
              CBUFFER = "tetra4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (BULK_STRAIN(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MM3 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MM3   
              CBUFFER = "tria3"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (BULK_STRAIN(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MP3 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MP3   
              CBUFFER = "tria3"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (BULK_STRAIN(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MM4 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MM4   
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (BULK_STRAIN(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MP4 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MP4 
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (BULK_STRAIN(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MTR .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MTR   
              CBUFFER = "bar2"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (BULK_STRAIN(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

          ENDIF
        ENDDO
!!
!! Close out this-time-step block.
!!
        CBUFFER = "END TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')

      ENDIF
!!
!! Strain Energy Density...
!!
      IOERROR = .TRUE.
      OPEN
     &  (
     &  UNIT     =  IO_UNIT%LEGO,
     &  FILE     = 'fmaego.data.esed',
     &  STATUS   = 'UNKNOWN',
#ifdef _G95_
     &  FORM     = 'UNFORMATTED', ACCESS='STREAM',  !  G95
#endif
#ifdef _CVF_NT_
     &  FORM     = 'BINARY', CONVERT='BIG_ENDIAN',  !  CVF
#endif
#ifdef LANGUAGE_FORTRAN90
     &  FORM     = 'BINARY',                        !  Pathf95
#endif
     &  POSITION = 'APPEND',
     &  ERR      =  704
     &  )
      IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 704  IF (IOERROR) THEN
        CALL USER_MESSAGE
     &    (
     &    MSGL//'FATAL'//
     &    MSGL//'WRITE_ESG_DATA_FILE.007.04'//
     &    MSGL//'Unable To Execute OPEN On: '//'fmaego.data.esed'
     &    )
      ELSE
!!
!! Initialize sequential Material part counter.
!!
        MPart = 0
!!
!! Start this-time-step block.
!!
        CBUFFER = "BEGIN TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CBUFFER = "Element Strain Energy Density Results"
        WRITE (IO_UNIT%LEGO) CBUFFER
!!
!! Loop on Material parts.
!!
        DO M = 1,NUMMT
!!
!! The function ELEMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! elements-used array NELUSED and vtk output-index translation array NPTNOWI for 
!! this material. If material M is used, the function returns 1, otherwise 0. 
!!
          IF (ELEMENTS_AND_NODES_USED(M) .GT. 0) THEN

            NNodes = MATERIAL(M)%NNodes
            NElems = MATERIAL(M)%NElems

            CBUFFER = "part"
            WRITE (IO_UNIT%LEGO) CBUFFER
            MPart = MPart + 1
            WRITE (IO_UNIT%LEGO) MPart

            Nend = 0

            IF (MHX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MHX
              CBUFFER = "hexa8"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRAIN_ENERGY_DENSITY(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MPX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MPX
              CBUFFER = "penta6"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRAIN_ENERGY_DENSITY(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MPY .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MPY   
              CBUFFER = "pyramid5"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRAIN_ENERGY_DENSITY(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MTX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MTX   
              CBUFFER = "tetra4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRAIN_ENERGY_DENSITY(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MM3 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MM3   
              CBUFFER = "tria3"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRAIN_ENERGY_DENSITY(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MP3 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MP3   
              CBUFFER = "tria3"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRAIN_ENERGY_DENSITY(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MM4 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MM4   
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRAIN_ENERGY_DENSITY(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MP4 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MP4 
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRAIN_ENERGY_DENSITY(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MTR .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MTR   
              CBUFFER = "bar2"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (STRAIN_ENERGY_DENSITY(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

          ENDIF
        ENDDO
!!
!! Close out this-time-step block.
!!
        CBUFFER = "END TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')

      ENDIF
!!
!! Pressure...
!!
      IOERROR = .TRUE.
      OPEN
     &  (
     &  UNIT     =  IO_UNIT%LEGO,
     &  FILE     = 'fmaego.data.eprs',
     &  STATUS   = 'UNKNOWN',
#ifdef _G95_
     &  FORM     = 'UNFORMATTED', ACCESS='STREAM',  !  G95
#endif
#ifdef _CVF_NT_
     &  FORM     = 'BINARY', CONVERT='BIG_ENDIAN',  !  CVF
#endif
#ifdef LANGUAGE_FORTRAN90
     &  FORM     = 'BINARY',                        !  Pathf95
#endif
     &  POSITION = 'APPEND',
     &  ERR      =  705
     &  )
      IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 705  IF (IOERROR) THEN
        CALL USER_MESSAGE
     &    (
     &    MSGL//'FATAL'//
     &    MSGL//'WRITE_ESG_DATA_FILE.007.05'//
     &    MSGL//'Unable To Execute OPEN On: '//'fmaego.data.eprs'
     &    )
      ELSE
!!
!! Initialize sequential Material part counter.
!!
        MPart = 0
!!
!! Start this-time-step block.
!!
        CBUFFER = "BEGIN TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CBUFFER = "Element Pressure Density Results"
        WRITE (IO_UNIT%LEGO) CBUFFER
!!
!! Loop on Material parts.
!!
        DO M = 1,NUMMT
!!
!! The function ELEMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! elements-used array NELUSED and vtk output-index translation array NPTNOWI for 
!! this material. If material M is used, the function returns 1, otherwise 0. 
!!
          IF (ELEMENTS_AND_NODES_USED(M) .GT. 0) THEN

            NNodes = MATERIAL(M)%NNodes
            NElems = MATERIAL(M)%NElems

            CBUFFER = "part"
            WRITE (IO_UNIT%LEGO) CBUFFER
            MPart = MPart + 1
            WRITE (IO_UNIT%LEGO) MPart

            Nend = 0

            IF (MHX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MHX
              CBUFFER = "hexa8"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (PRESSURE(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MPX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MPX
              CBUFFER = "penta6"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (PRESSURE(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MPY .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MPY   
              CBUFFER = "pyramid5"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (PRESSURE(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MTX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MTX   
              CBUFFER = "tetra4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (PRESSURE(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MM3 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MM3   
              CBUFFER = "tria3"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (PRESSURE(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MP3 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MP3   
              CBUFFER = "tria3"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (PRESSURE(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MM4 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MM4   
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (PRESSURE(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MP4 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MP4 
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (PRESSURE(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MTR .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MTR   
              CBUFFER = "bar2"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (PRESSURE(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

          ENDIF
        ENDDO
!!
!! Close out this-time-step block.
!!
        CBUFFER = "END TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')

      ENDIF
!!
!! Effective stress (deviatoric stress magnitude)...
!!
      IOERROR = .TRUE.
      OPEN
     &  (
     &  UNIT     =  IO_UNIT%LEGO,
     &  FILE     = 'fmaego.data.edev',
     &  STATUS   = 'UNKNOWN',
#ifdef _G95_
     &  FORM     = 'UNFORMATTED', ACCESS='STREAM',  !  G95
#endif
#ifdef _CVF_NT_
     &  FORM     = 'BINARY', CONVERT='BIG_ENDIAN',  !  CVF
#endif
#ifdef LANGUAGE_FORTRAN90
     &  FORM     = 'BINARY',                        !  Pathf95
#endif
     &  POSITION = 'APPEND',
     &  ERR      =  706
     &  )
      IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 706  IF (IOERROR) THEN
        CALL USER_MESSAGE
     &    (
     &    MSGL//'FATAL'//
     &    MSGL//'WRITE_ESG_DATA_FILE.007.06'//
     &    MSGL//'Unable To Execute OPEN On: '//'fmaego.data.edev'
     &    )
      ELSE
!!
!! Initialize sequential Material part counter.
!!
        MPart = 0
!!
!! Start this-time-step block.
!!
        CBUFFER = "BEGIN TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CBUFFER = "Element Effective Stress Results"
        WRITE (IO_UNIT%LEGO) CBUFFER
!!
!! Loop on Material parts.
!!
        DO M = 1,NUMMT
!!
!! The function ELEMENTS_AND_NODES_USED(*) gen's nodal-points-used array NPTUSED, 
!! elements-used array NELUSED and vtk output-index translation array NPTNOWI for 
!! this material. If material M is used, the function returns 1, otherwise 0. 
!!
          IF (ELEMENTS_AND_NODES_USED(M) .GT. 0) THEN

            NNodes = MATERIAL(M)%NNodes
            NElems = MATERIAL(M)%NElems

            CBUFFER = "part"
            WRITE (IO_UNIT%LEGO) CBUFFER
            MPart = MPart + 1
            WRITE (IO_UNIT%LEGO) MPart

            Nend = 0

            IF (MHX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MHX
              CBUFFER = "hexa8"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (EFFECTIVE_STRESS(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MPX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MPX
              CBUFFER = "penta6"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (EFFECTIVE_STRESS(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MPY .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MPY   
              CBUFFER = "pyramid5"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (EFFECTIVE_STRESS(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MTX .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MTX   
              CBUFFER = "tetra4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (EFFECTIVE_STRESS(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MM3 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MM3   
              CBUFFER = "tria3"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (EFFECTIVE_STRESS(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MP3 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MP3   
              CBUFFER = "tria3"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (EFFECTIVE_STRESS(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MM4 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MM4   
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (EFFECTIVE_STRESS(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MP4 .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MP4 
              CBUFFER = "quad4"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (EFFECTIVE_STRESS(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

            IF (MTR .GT. 0) THEN
              Nbgn = Nend + 1
              Nend = Nend + MTR   
              CBUFFER = "bar2"
              WRITE (IO_UNIT%LEGO) CBUFFER
              WRITE (IO_UNIT%LEGO) (EFFECTIVE_STRESS(NELUSED(n)), n = Nbgn,Nend)
            ENDIF

          ENDIF
        ENDDO
!!
!! Close out this-time-step block.
!!
        CBUFFER = "END TIME STEP"
        WRITE (IO_UNIT%LEGO) CBUFFER

        CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')

      ENDIF
!!
!! =========FOURTH============================================================
!! Write a *.case file to allow EnSight and ParaView to access the nodal point
!! and element (cell) results.
!! 
!! Open an ESG case-file.
!!
      IOERROR = .TRUE.
      OPEN
     &  (
     &  UNIT    =  IO_UNIT%LEGO,
     &  FILE    = 'fmaego.data.case',
     &  STATUS  = 'UNKNOWN',
     &  FORM    = 'FORMATTED',
     &  ERR     =  800
     &  )
      IOERROR = .FALSE.
!!
!! Fatal error exit for failed OPEN operation.
!!
 800  IF (IOERROR) THEN
        CALL USER_MESSAGE
     &    (
     &    MSGL//'FATAL'//
     &    MSGL//'WRITE_ESG_DATA_FILE.010.03'//
     &    MSGL//'Unable To Execute OPEN On: '//'fmaego.data.case'
     &    )
      ELSE

        WRITE (IO_UNIT%LEGO,'(A)') "FORMAT"
        WRITE (IO_UNIT%LEGO,'(A)') "type: ensight gold"

        !                                                   ts fs
        WRITE (IO_UNIT%LEGO,'(A)') "GEOMETRY"
        WRITE (IO_UNIT%LEGO,'(A)') "model:                   1  1  fmaego.data.geom    change_coords_only"

        !                                                   ts fs
        WRITE (IO_UNIT%LEGO,'(A)') "VARIABLE"
        WRITE (IO_UNIT%LEGO,'(A)') "vector per node:         1  1  Displacement          fmaego.data.ndis"
        WRITE (IO_UNIT%LEGO,'(A)') "vector per node:         1  1  Velocity              fmaego.data.nvel"
        WRITE (IO_UNIT%LEGO,'(A)') "vector per node:         1  1  Acceleration          fmaego.data.nacc"
        WRITE (IO_UNIT%LEGO,'(A)') "scalar per element:      1  1  Material              fmaego.data.emat"
        WRITE (IO_UNIT%LEGO,'(A)') "scalar per element:      1  1  Pressure              fmaego.data.eprs"
        WRITE (IO_UNIT%LEGO,'(A)') "scalar per element:      1  1  Bulk-Strain           fmaego.data.elnv"
        WRITE (IO_UNIT%LEGO,'(A)') "scalar per element:      1  1  Strain-Energy-Density fmaego.data.esed"
        WRITE (IO_UNIT%LEGO,'(A)') "scalar per element:      1  1  Effective-Stress      fmaego.data.edev"
        WRITE (IO_UNIT%LEGO,'(A)') "tensor symm per element: 1  1  Stress                fmaego.data.estr"

        !                                                   ts
        WRITE (IO_UNIT%LEGO,'(A)') "TIME"
        WRITE (IO_UNIT%LEGO,'(A)') "time set:                1"
        WRITE (IO_UNIT%LEGO,'(A)') "number of steps:     "//NEGO
        WRITE (IO_UNIT%LEGO,'(A,12X,2(1PE24.12)/(3(1PE24.12)))') "time values:",(TIME_VALUES(n), n= 1,TStep)

        !                                                      fs
        WRITE (IO_UNIT%LEGO,'(A)') "FILE"
        WRITE (IO_UNIT%LEGO,'(A)') "file set:                   1"
        WRITE (IO_UNIT%LEGO,'(A)') "number of steps:        "//NEGO

        CLOSE (UNIT=IO_UNIT%LEGO, STATUS='KEEP')
      ENDIF

      RETURN

      CONTAINS
!!
      FUNCTION ELEMENTS_AND_NODES_USED ( M )
!!
!! Copyright (c) by FMA Development, LLC, 26-MAY-2004
!!
!! Purpose: Mark elements and nodes used, count number of unique elements 
!! and nodes used by the elements, and create index translation arrays 
!! for material "M".
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: M  !  Material Index
!!
!! Function return value.
      INTEGER(4) :: ELEMENTS_AND_NODES_USED
!!
!! Local Variables.
      INTEGER(ITYPE), PARAMETER :: IUNITY(8) = (/1,1,1,1,1,1,1,1/)
      INTEGER(ITYPE) :: K,N
!!
!! Clear marker/index-sequence and translation arrays.
!!
      NELUSED = 0
      NPTUSED = 0
      NPTNOWI = 0
!!
!! Initialize the program's element-type counters
!!
      MHX = 0;  MPX = 0;  MPY = 0;  MTX = 0;  
      MM3 = 0;  MP3 = 0;  MM4 = 0;  MP4 = 0;  
      MTR = 0;
!!
!! Scan all elements for the elements using MATERIAL(M).
!! Note: The element scan order here must match that in
!! the other extraction routines: NUMHX, NUMPX, NUMPY,
!! NUMTX, NUMM3, NUMP3, NUMM4, NUMP4, NUMTR.
!!
      K = 0
      ELEMENTS_AND_NODES_USED = 0

      DO N = 1,NUMHX
        K = K + 1
        IF (HEXAH(N)%PAR%MatID .EQ. M) THEN
          MHX = MHX + 1
          NELUSED(K) = 1
          ELEMENTS_AND_NODES_USED = 1
          NPTUSED(HEXAH(N)%PAR%IX(1:8)) = IUNITY(1:8)
        ENDIF
      ENDDO

      DO N = 1,NUMPX
        K = K + 1
        IF (PENTA(N)%PAR%MatID .EQ. M) THEN
          MPX = MPX + 1
          NELUSED(K) = 1
          ELEMENTS_AND_NODES_USED = 1
          NPTUSED(PENTA(N)%PAR%IX(1:6)) = IUNITY(1:6)
        ENDIF
      ENDDO

      DO N = 1,NUMPY
        K = K + 1
        IF (PYRAMID(N)%PAR%MatID .EQ. M) THEN
          MPY = MPY + 1
          NELUSED(K) = 1
          ELEMENTS_AND_NODES_USED = 1
          NPTUSED(PYRAMID(N)%PAR%IX(1:5)) = IUNITY(1:5)
        ENDIF
      ENDDO

      DO N = 1,NUMTX
        K = K + 1
        IF (TETRA(N)%PAR%MatID .EQ. M) THEN
          MTX = MTX + 1
          NELUSED(K) = 1
          ELEMENTS_AND_NODES_USED = 1

! ParaView does not know how to handle 8-node tetrahedrons.
!          Iform = SECTION_3D( TETRA(N)%PAR%SecID )%Iform
!          IF (Iform .EQ. Eight_Nodes) Then
!            NPTUSED(TETRA(N)%PAR%IX(1:8)) = IUNITY(1:8)
!          ELSE
!            NPTUSED(TETRA(N)%PAR%IX(1:4)) = IUNITY(1:4)
!          ENDIF

          NPTUSED(TETRA(N)%PAR%IX(1:4)) = IUNITY(1:4)
        ENDIF
      ENDDO

      DO N = 1,NUMM3
        K = K + 1
        IF (MEMBT(N)%PAR%MatID .EQ. M) THEN
          MM3 = MM3 + 1
          NELUSED(K) = 1
          ELEMENTS_AND_NODES_USED = 1
          NPTUSED(MEMBT(N)%PAR%IX(1:3)) = IUNITY(1:3)
        ENDIF
      ENDDO

      DO N = 1,NUMP3
        K = K + 1
        IF (PLATT(N)%PAR%MatID .EQ. M) THEN
          MP3 = MP3 + 1
          NELUSED(K) = 1
          ELEMENTS_AND_NODES_USED = 1
          NPTUSED(PLATT(N)%PAR%IX(1:3)) = IUNITY(1:3)
        ENDIF
      ENDDO

      DO N = 1,NUMM4
        K = K + 1
        IF (MEMBQ(N)%PAR%MatID .EQ. M) THEN
          MM4 = MM4 + 1
          NELUSED(K) = 1
          ELEMENTS_AND_NODES_USED = 1
          NPTUSED(MEMBQ(N)%PAR%IX(1:4)) = IUNITY(1:4)
        ENDIF
      ENDDO

      DO N = 1,NUMP4
        K = K + 1
        IF (PLATQ(N)%PAR%MatID .EQ. M) THEN
          MP4 = MP4 + 1
          NELUSED(K) = 1
          ELEMENTS_AND_NODES_USED = 1
          NPTUSED(PLATQ(N)%PAR%IX(1:4)) = IUNITY(1:4)
        ENDIF
      ENDDO

      DO N = 1,NUMTR
        K = K + 1
        IF (TRUSS(N)%PAR%MatID .EQ. M) THEN
          MTR = MTR + 1
          NELUSED(K) = 1
          ELEMENTS_AND_NODES_USED = 1
          NPTUSED(TRUSS(N)%PAR%IX(1:2)) = IUNITY(1:2)
        ENDIF
      ENDDO
!!
!! Get count of elements and nodes used by this material.
!!
      MATERIAL(M)%NElems = SUM ( NELUSED(1:NUMXL) )
      MATERIAL(M)%NNodes = SUM ( NPTUSED(1:NUMNP) )

!!    write (IO_UNIT%LELO,*) ' *** Material Number:',M
!!    write (IO_UNIT%LELO,*) '              NElems:',MATERIAL(M)%NElems
!!    write (IO_UNIT%LELO,*) '              NNodes:',MATERIAL(M)%NNodes
!!    write (IO_UNIT%LELO,*) '    NELUSED(1:NUMXL):',(NELUSED(n), n=1,NUMXL)
!!    write (IO_UNIT%LELO,*) '    NPTUSED(1:NUMNP):',(NPTUSED(n), n=1,NUMNP)

!!
!! Create nodal point index map NPTNOWI for use with 
!! element connectivity n-tuples. (Maps program index 
!! to index written to esg files for this material.)
!!
!! Also, covert NPTUSED into a sequential index map.
!!
      K = 0
      DO N = 1,NUMNP
        IF (NPTUSED(N) .EQ. 1) THEN
          K = K + 1
          NPTNOWI(N) = K
          NPTUSED(K) = N
        ENDIF
      ENDDO
!!
!! Covert NELUSED into a sequential index map.
!!
      K= 0 
      DO N = 1,NUMXL
        IF (NELUSED(N) .EQ. 1) THEN
          K = K + 1
          NELUSED(K) = N
        ENDIF
      ENDDO

      RETURN
      END FUNCTION ELEMENTS_AND_NODES_USED
!!
      FUNCTION ELEMENT_ESG_TYPE( NEL )
!!
!! Copyright (c) by FMA Development, LLC, 26-MAY-2004
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: NEL   ! Nth element
!!
!! Function return value.
      CHARACTER(8) :: ELEMENT_ESG_TYPE
!!
!! Local variables.
      LOGICAL,        SAVE :: FIRST = .TRUE.
      INTEGER(ITYPE), SAVE :: NHX,NPX,NPY,NTX,NM3,NP3,NM4,NP4,NTR,N

      CHARACTER(8), SAVE :: ESG_VERTEX      = "point"
      CHARACTER(8), SAVE :: ESG_LINE        = "bar2"
      CHARACTER(8), SAVE :: ESG_TRIANGLE    = "tria3"
      CHARACTER(8), SAVE :: ESG_QUAD        = "quad4"
      CHARACTER(8), SAVE :: ESG_TETRAHEDRON = "tetra4"
      CHARACTER(8), SAVE :: ESG_HEXAHEDRON  = "hexa8"
      CHARACTER(8), SAVE :: ESG_WEDGE       = "penta6"
      CHARACTER(8), SAVE :: ESG_PYRAMID     = "pyramid5"
!!
      IF (FIRST) THEN
        NHX = NUMHX
        NPX = NHX + NUMPX
        NPY = NPX + NUMPY
        NTX = NPY + NUMTX
        NM3 = NTX + NUMM3
        NP3 = NM3 + NUMP3
        NM4 = NP3 + NUMM4
        NP4 = NM4 + NUMP4
        NTR = NP4 + NUMTR
        FIRST = .FALSE.
      ENDIF
!!
      ELEMENT_ESG_TYPE = " "
      IF (NEL .LE. NHX) THEN
        N = NEL
        ELEMENT_ESG_TYPE = ESG_HEXAHEDRON  ! HEXAH(N)%PAR%EleID
      ELSE IF (NEL .LE. NPX) THEN
        N = NEL - NHX
        ELEMENT_ESG_TYPE = ESG_WEDGE       ! PENTA(N)%PAR%EleID
      ELSE IF (NEL .LE. NPY) THEN
        N = NEL - NPX
        ELEMENT_ESG_TYPE = ESG_PYRAMID     ! PYRAMID(N)%PAR%EleID
      ELSE IF (NEL .LE. NTX) THEN
        N = NEL - NPY
        ELEMENT_ESG_TYPE = ESG_TETRAHEDRON ! TETRA(N)%PAR%EleID
      ELSE IF (NEL .LE. NM3) THEN
        N = NEL - NTX
        ELEMENT_ESG_TYPE = ESG_TRIANGLE    ! MEMBT(N)%PAR%EleID
      ELSE IF (NEL .LE. NP3) THEN
        N = NEL - NM3
        ELEMENT_ESG_TYPE = ESG_TRIANGLE    ! PLATT(N)%PAR%EleID
      ELSE IF (NEL .LE. NM4) THEN
        N = NEL - NP3
        ELEMENT_ESG_TYPE = ESG_QUAD        ! MEMBQ(N)%PAR%EleID
      ELSE IF (NEL .LE. NP4) THEN
        N = NEL - NM4
        ELEMENT_ESG_TYPE = ESG_QUAD        ! PLATQ(N)%PAR%EleID
      ELSE IF (NEL .LE. NTR) THEN
        N = NEL - NP4
        ELEMENT_ESG_TYPE = ESG_LINE        ! TRUSS(N)%PAR%EleID
      ENDIF

      RETURN
      END FUNCTION ELEMENT_ESG_TYPE
!!
      FUNCTION EleID_DATA ( NEL )
!!
!! Copyright (c) by FMA Development, LLC, 12-APR-1991 18:59:44 
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: NEL   ! Nth element
!!
!! Function return value.
      INTEGER(4) :: EleID_DATA
!!
!! Local variables.
      LOGICAL,        SAVE :: FIRST = .TRUE.
      INTEGER(ITYPE), SAVE :: NHX,NPX,NPY,NTX,NM3,NP3,NM4,NP4,NTR,N
!!
      IF (FIRST) THEN
        NHX = NUMHX
        NPX = NHX + NUMPX
        NPY = NPX + NUMPY
        NTX = NPY + NUMTX
        NM3 = NTX + NUMM3
        NP3 = NM3 + NUMP3
        NM4 = NP3 + NUMM4
        NP4 = NM4 + NUMP4
        NTR = NP4 + NUMTR
        FIRST = .FALSE.
      ENDIF
!!
      EleID_DATA = 0
      IF (NEL .LE. NHX) THEN
        N = NEL
        EleID_DATA = HEXAH(N)%PAR%EleID
      ELSE IF (NEL .LE. NPX) THEN
        N = NEL - NHX
        EleID_DATA = PENTA(N)%PAR%EleID
      ELSE IF (NEL .LE. NPY) THEN
        N = NEL - NPX
        EleID_DATA = PYRAMID(N)%PAR%EleID
      ELSE IF (NEL .LE. NTX) THEN
        N = NEL - NPY
        EleID_DATA = TETRA(N)%PAR%EleID
      ELSE IF (NEL .LE. NM3) THEN
        N = NEL - NTX
        EleID_DATA = MEMBT(N)%PAR%EleID
      ELSE IF (NEL .LE. NP3) THEN
        N = NEL - NM3
        EleID_DATA = PLATT(N)%PAR%EleID
      ELSE IF (NEL .LE. NM4) THEN
        N = NEL - NP3
        EleID_DATA = MEMBQ(N)%PAR%EleID
      ELSE IF (NEL .LE. NP4) THEN
        N = NEL - NM4
        EleID_DATA = PLATQ(N)%PAR%EleID
      ELSE IF (NEL .LE. NTR) THEN
        N = NEL - NP4
        EleID_DATA = TRUSS(N)%PAR%EleID
      ENDIF

      RETURN
      END FUNCTION EleID_DATA
!!
      FUNCTION ELEMENT_N_TUPLE( NEL )
!!
!! Copyright (c) by FMA Development, LLC, 26-MAY-2004
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: NEL   ! Nth element
!!
!! Function return value.
      INTEGER(ITYPE) :: ELEMENT_N_TUPLE
!!
!! Local variables.
      LOGICAL,        SAVE :: FIRST = .TRUE.
      INTEGER(ITYPE), SAVE :: NHX,NPX,NPY,NTX,NM3,NP3,NM4,NP4,NTR,N
!!
      IF (FIRST) THEN
        NHX = NUMHX
        NPX = NHX + NUMPX
        NPY = NPX + NUMPY
        NTX = NPY + NUMTX
        NM3 = NTX + NUMM3
        NP3 = NM3 + NUMP3
        NM4 = NP3 + NUMM4
        NP4 = NM4 + NUMP4
        NTR = NP4 + NUMTR
        FIRST = .FALSE.
      ENDIF
!!
      ELEMENT_N_TUPLE = 0
      IF (NEL .LE. NHX) THEN
        N = NEL
        ELEMENT_N_TUPLE = 8
        NTUPLE(1:8) = NPTNOWI(HEXAH(N)%PAR%IX(1:8))
      ELSE IF (NEL .LE. NPX) THEN
        N = NEL - NHX
        ELEMENT_N_TUPLE = 6
        NTUPLE(1:6) = NPTNOWI(PENTA(N)%PAR%IX(1:6))
      ELSE IF (NEL .LE. NPY) THEN
        N = NEL - NPX
        ELEMENT_N_TUPLE = 5
        NTUPLE(1:5) = NPTNOWI(PYRAMID(N)%PAR%IX(1:5))
      ELSE IF (NEL .LE. NTX) THEN
        N = NEL - NPY
        ELEMENT_N_TUPLE = 4
        NTUPLE(1:4) = NPTNOWI(TETRA(N)%PAR%IX(1:4))
      ELSE IF (NEL .LE. NM3) THEN
        N = NEL - NTX
        ELEMENT_N_TUPLE = 3
        NTUPLE(1:3) = NPTNOWI(MEMBT(N)%PAR%IX(1:3))
      ELSE IF (NEL .LE. NP3) THEN
        N = NEL - NM3
        ELEMENT_N_TUPLE = 3
        NTUPLE(1:3) = NPTNOWI(PLATT(N)%PAR%IX(1:3))
      ELSE IF (NEL .LE. NM4) THEN
        N = NEL - NP3
        ELEMENT_N_TUPLE = 4
        NTUPLE(1:4) = NPTNOWI(MEMBQ(N)%PAR%IX(1:4))
      ELSE IF (NEL .LE. NP4) THEN
        N = NEL - NM4
        ELEMENT_N_TUPLE = 4
        NTUPLE(1:4) = NPTNOWI(PLATQ(N)%PAR%IX(1:4))
      ELSE IF (NEL .LE. NTR) THEN
        N = NEL - NP4
        ELEMENT_N_TUPLE = 2
        NTUPLE(1:2) = NPTNOWI(TRUSS(N)%PAR%IX(1:2))
      ENDIF

      RETURN
      END FUNCTION ELEMENT_N_TUPLE
!!
      FUNCTION MatID_DATA ( NEL )
!!
!! Copyright (c) by FMA Development, LLC, 12-APR-1991 18:59:44 
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: NEL   ! Nth element
!!
!! Function return value.
      INTEGER(ITYPE) :: MatID_DATA
!!
!! Local variables.
      LOGICAL,        SAVE :: FIRST = .TRUE.
      INTEGER(ITYPE), SAVE :: NHX,NPX,NPY,NTX,NM3,NP3,NM4,NP4,NTR,N
!!
      IF (FIRST) THEN
        NHX = NUMHX
        NPX = NHX + NUMPX
        NPY = NPX + NUMPY
        NTX = NPY + NUMTX
        NM3 = NTX + NUMM3
        NP3 = NM3 + NUMP3
        NM4 = NP3 + NUMM4
        NP4 = NM4 + NUMP4
        NTR = NP4 + NUMTR
        FIRST = .FALSE.
      ENDIF
!!
      MatID_DATA = 0
      IF (NEL .LE. NHX) THEN
        N = NEL
        MatID_DATA = HEXAH(N)%PAR%MatID
      ELSE IF (NEL .LE. NPX) THEN
        N = NEL - NHX
        MatID_DATA = PENTA(N)%PAR%MatID
      ELSE IF (NEL .LE. NPY) THEN
        N = NEL - NPX
        MatID_DATA = PYRAMID(N)%PAR%MatID
      ELSE IF (NEL .LE. NTX) THEN
        N = NEL - NPY
        MatID_DATA = TETRA(N)%PAR%MatID
      ELSE IF (NEL .LE. NM3) THEN
        N = NEL - NTX
        MatID_DATA = MEMBT(N)%PAR%MatID
      ELSE IF (NEL .LE. NP3) THEN
        N = NEL - NM3
        MatID_DATA = PLATT(N)%PAR%MatID
      ELSE IF (NEL .LE. NM4) THEN
        N = NEL - NP3
        MatID_DATA = MEMBQ(N)%PAR%MatID
      ELSE IF (NEL .LE. NP4) THEN
        N = NEL - NM4
        MatID_DATA = PLATQ(N)%PAR%MatID
      ELSE IF (NEL .LE. NTR) THEN
        N = NEL - NP4
        MatID_DATA = TRUSS(N)%PAR%MatID
      ENDIF

      RETURN
      END FUNCTION MatID_DATA
!!
      FUNCTION STRESS_DATA ( NEL, I )
!!
!! Copyright (c) by FMA Development, LLC, 8-APR-1991 20:33:30 
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: NEL  ! Nth element
      INTEGER(ITYPE), INTENT(IN) :: I    ! Stress component requested
!!
!! Function return value.
      REAL :: STRESS_DATA
!!
!! Local variables.
      LOGICAL,        SAVE :: FIRST = .TRUE.
      INTEGER(ITYPE), SAVE :: NHX,NPX,NPY,NTX,NM3,NP3,NM4,NP4,NTR,N
      INTEGER(ITYPE)       :: NDEX(4)
!!
      IF (FIRST) THEN
        NHX = NUMHX
        NPX = NHX + NUMPX
        NPY = NPX + NUMPY
        NTX = NPY + NUMTX
        NM3 = NTX + NUMM3
        NP3 = NM3 + NUMP3
        NM4 = NP3 + NUMM4
        NP4 = NM4 + NUMP4
        NTR = NP4 + NUMTR
        FIRST = .FALSE.
      ENDIF
!!
      STRESS_DATA = 0.0
      IF (NEL .LE. NHX) THEN
        N = NEL
        STRESS_DATA = HEXAH(N)%RES%Stress(I)
      ELSE IF (NEL .LE. NPX) THEN
        N = NEL - NHX
        STRESS_DATA = PENTA(N)%RES%Stress(I)
      ELSE IF (NEL .LE. NPY) THEN
        N = NEL - NPX
        STRESS_DATA = PYRAMID(N)%RES%Stress(I)
      ELSE IF (NEL .LE. NTX) THEN
        N = NEL - NPY
        Iform = SECTION_3D( TETRA(N)%PAR%SecID )%Iform
        IF (Iform .EQ. Nodal_Based) Then
          NDEX = TETRA(N)%RES%IN
          STRESS_DATA = P4TH * SUM( TETNP(NDEX)%Stress(I) )
        ELSE
          STRESS_DATA = TETRA(N)%RES%Stress(I)
        ENDIF
      ELSE IF (NEL .LE. NM3) THEN
        N = NEL - NTX
        IF (I .LE. 2) THEN
          STRESS_DATA = MEMBT(N)%RES%Stress(I)
        ELSE IF (I .EQ. 4) THEN
          STRESS_DATA = MEMBT(N)%RES%Stress(3)
        ENDIF
      ELSE IF (NEL .LE. NP3) THEN
        N = NEL - NM3
        Ist = PLATT(N)%PAR%Ist
        STRESS_DATA = STRESS(I,Ist)
      ELSE IF (NEL .LE. NM4) THEN
        N = NEL - NP3
        IF (I .LE. 2) THEN
          STRESS_DATA = MEMBQ(N)%RES%Stress(I)
        ELSE IF (I .EQ. 4) THEN
          STRESS_DATA = MEMBQ(N)%RES%Stress(3)
        ENDIF
      ELSE IF (NEL .LE. NP4) THEN
        N = NEL - NM4
        Ist = PLATQ(N)%PAR%Ist
        STRESS_DATA = STRESS(I,Ist)
      ELSE IF (NEL .LE. NTR) THEN
        N = NEL - NP4
        IF (I .EQ. 1) THEN
          STRESS_DATA = TRUSS(N)%RES%Stress
        ELSE
          STRESS_DATA = 0.0
        ENDIF
      ENDIF

      RETURN
      END FUNCTION STRESS_DATA
!!
      FUNCTION BULK_STRAIN ( NEL )
!!
!! Copyright (c) by FMA Development, LLC, 12-APR-1991 18:59:44 
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: NEL  ! Nth element
!!
!! Function return value.
      REAL :: BULK_STRAIN
!!
!! Local variables.
      LOGICAL,        SAVE :: FIRST = .TRUE.
      INTEGER(ITYPE), SAVE :: NHX,NPX,NPY,NTX,NM3,NP3,NM4,NP4,NTR,N
!!
!! Function reference.
      REAL(RTYPE), EXTERNAL :: LOG_POS
!!
      IF (FIRST) THEN
        NHX = NUMHX
        NPX = NHX + NUMPX
        NPY = NPX + NUMPY
        NTX = NPY + NUMTX
        NM3 = NTX + NUMM3
        NP3 = NM3 + NUMP3
        NM4 = NP3 + NUMM4
        NP4 = NM4 + NUMP4
        NTR = NP4 + NUMTR
        FIRST = .FALSE.
      ENDIF
!!
      BULK_STRAIN = 0.0
      IF (NEL .LE. NHX) THEN
        N = NEL
        BULK_STRAIN = LOG_POS(HEXAH(N)%RES%Volume/HEXAH(N)%PAR%Volume)
      ELSE IF (NEL .LE. NPX) THEN
        N = NEL - NHX
        BULK_STRAIN = LOG_POS(PENTA(N)%RES%Volume/PENTA(N)%PAR%Volume)
      ELSE IF (NEL .LE. NPY) THEN
        N = NEL - NPX
        BULK_STRAIN = LOG_POS(PYRAMID(N)%RES%Volume/PYRAMID(N)%PAR%Volume)
      ELSE IF (NEL .LE. NTX) THEN
        N = NEL - NPY
        BULK_STRAIN = LOG_POS(TETRA(N)%RES%Volume/TETRA(N)%PAR%Volume)
      ELSE IF (NEL .LE. NM3) THEN
        N = NEL - NTX
        BULK_STRAIN = LOG_POS(MEMBT(N)%RES%Area/MEMBT(N)%PAR%Area)
      ELSE IF (NEL .LE. NP3) THEN
        N = NEL - NM3
        BULK_STRAIN = LOG_POS(PLATT(N)%RES%Area/PLATT(N)%PAR%Area)
      ELSE IF (NEL .LE. NM4) THEN
        N = NEL - NP3
        BULK_STRAIN = LOG_POS(MEMBQ(N)%RES%Area/MEMBQ(N)%PAR%Area)
      ELSE IF (NEL .LE. NP4) THEN
        N = NEL - NM4
        BULK_STRAIN = LOG_POS(PLATQ(N)%RES%Area/PLATQ(N)%PAR%Area)
      ELSE IF (NEL .LE. NTR) THEN
        N = NEL - NP4
        BULK_STRAIN = LOG_POS(TRUSS(N)%RES%Length/TRUSS(N)%PAR%Length)
      ENDIF

      RETURN
      END FUNCTION BULK_STRAIN
!!
      FUNCTION STRAIN_ENERGY_DENSITY ( NEL )
!!
!! Copyright (c) by FMA Development, LLC, 12-APR-1991 18:59:44 
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: NEL   ! Nth element
!!
!! Function return value.
      REAL :: STRAIN_ENERGY_DENSITY
!!
!! Local variables.
      LOGICAL,        SAVE :: FIRST = .TRUE.
      INTEGER(ITYPE), SAVE :: NHX,NPX,NPY,NTX,NM3,NP3,NM4,NP4,NTR,N
      INTEGER(ITYPE)       :: NDEX(4)
!!
      IF (FIRST) THEN
        NHX = NUMHX
        NPX = NHX + NUMPX
        NPY = NPX + NUMPY
        NTX = NPY + NUMTX
        NM3 = NTX + NUMM3
        NP3 = NM3 + NUMP3
        NM4 = NP3 + NUMM4
        NP4 = NM4 + NUMP4
        NTR = NP4 + NUMTR
        FIRST = .FALSE.
      ENDIF
!!
      STRAIN_ENERGY_DENSITY = 0.0
      IF (NEL .LE. NHX) THEN
        N = NEL
        STRAIN_ENERGY_DENSITY = HEXAH(N)%RES%Str_Eng
      ELSE IF (NEL .LE. NPX) THEN
        N = NEL - NHX
        STRAIN_ENERGY_DENSITY = PENTA(N)%RES%Str_Eng
      ELSE IF (NEL .LE. NPY) THEN
        N = NEL - NPX
        STRAIN_ENERGY_DENSITY = PYRAMID(N)%RES%Str_Eng
      ELSE IF (NEL .LE. NTX) THEN
        N = NEL - NPY
        Iform = SECTION_3D( TETRA(N)%PAR%SecID )%Iform
        IF (Iform .EQ. Nodal_Based) Then
          NDEX = TETRA(N)%RES%IN
          STRAIN_ENERGY_DENSITY = P4TH * SUM( TETNP(NDEX)%Str_Eng )
        ELSE
          STRAIN_ENERGY_DENSITY = TETRA(N)%RES%Str_Eng
        ENDIF
      ELSE IF (NEL .LE. NM3) THEN
        N = NEL - NTX
        STRAIN_ENERGY_DENSITY = MEMBT(N)%RES%Str_Eng
      ELSE IF (NEL .LE. NP3) THEN
        N = NEL - NM3
        STRAIN_ENERGY_DENSITY = PLATT(N)%RES%Str_Eng
      ELSE IF (NEL .LE. NM4) THEN
        N = NEL - NP3
        STRAIN_ENERGY_DENSITY = MEMBQ(N)%RES%Str_Eng
      ELSE IF (NEL .LE. NP4) THEN
        N = NEL - NM4
        STRAIN_ENERGY_DENSITY = PLATQ(N)%RES%Str_Eng
      ELSE IF (NEL .LE. NTR) THEN
        N = NEL - NP4
        STRAIN_ENERGY_DENSITY = TRUSS(N)%RES%Str_Eng
      ENDIF
!!
      RETURN
      END FUNCTION STRAIN_ENERGY_DENSITY
!!
      FUNCTION PRESSURE ( NEL )
!!
!! Copyright (c) by FMA Development, LLC, 8-APR-1991 20:33:30 
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: NEL   ! Nth element
!!
!! Function return value.
      REAL :: PRESSURE
!!
!! Local variables.
      LOGICAL,        SAVE :: FIRST = .TRUE.
      INTEGER(ITYPE), SAVE :: NHX,NPX,NPY,NTX,NM3,NP3,NM4,NP4,NTR,N
      INTEGER(ITYPE)       :: NDEX(4)
      REAL(RTYPE)          :: TEMP(4)
!!
      IF (FIRST) THEN
        NHX = NUMHX
        NPX = NHX + NUMPX
        NPY = NPX + NUMPY
        NTX = NPY + NUMTX
        NM3 = NTX + NUMM3
        NP3 = NM3 + NUMP3
        NM4 = NP3 + NUMM4
        NP4 = NM4 + NUMP4
        NTR = NP4 + NUMTR
        FIRST = .FALSE.
      ENDIF
!!
      PRESSURE = 0.0
      IF (NEL .LE. NHX) THEN
        N = NEL
        PRESSURE = P3RD * SUM( HEXAH(N)%RES%Stress(1:3) )
      ELSE IF (NEL .LE. NPX) THEN
        N = NEL - NHX
        PRESSURE = P3RD * SUM( PENTA(N)%RES%Stress(1:3) )
      ELSE IF (NEL .LE. NPY) THEN
        N = NEL - NPX
        PRESSURE = P3RD * SUM( PYRAMID(N)%RES%Stress(1:3) )
      ELSE IF (NEL .LE. NTX) THEN
        N = NEL - NPY
        Iform = SECTION_3D( TETRA(N)%PAR%SecID )%Iform
        IF (Iform .EQ. Nodal_Based) Then
          NDEX = TETRA(N)%RES%IN
          TEMP = (/ZERO,ZERO,ZERO,ZERO/)
          TEMP = TEMP + TETNP(NDEX)%Stress(1)
          TEMP = TEMP + TETNP(NDEX)%Stress(2)
          TEMP = TEMP + TETNP(NDEX)%Stress(3)
          PRESSURE = P12TH * SUM( TEMP ) 
        ELSE
          PRESSURE = P3RD * SUM( TETRA(N)%RES%Stress(1:3) )
        ENDIF
      ELSE IF (NEL .LE. NM3) THEN
        N = NEL - NTX
        PRESSURE = P3RD * SUM( MEMBT(N)%RES%Stress(1:2) )
      ELSE IF (NEL .LE. NP3) THEN
        N = NEL - NM3
        Ist = PLATT(N)%PAR%Ist
        PRESSURE = P3RD * SUM( STRESS(1:2,Ist) )
      ELSE IF (NEL .LE. NM4) THEN
        N = NEL - NP3
        PRESSURE = P3RD * SUM( MEMBQ(N)%RES%Stress(1:2) )
      ELSE IF (NEL .LE. NP4) THEN
        N = NEL - NM4
        Ist = PLATQ(N)%PAR%Ist
        PRESSURE = P3RD * SUM( STRESS(1:2,Ist) )
      ELSE IF (NEL .LE. NTR) THEN
        N = NEL - NP4
        PRESSURE = P3RD * TRUSS(N)%RES%Stress
      ENDIF

      RETURN
      END FUNCTION PRESSURE
!!
      FUNCTION EFFECTIVE_STRESS ( NEL )
!!
!! Copyright (c) by FMA Development, LLC, 8-APR-1991 20:33:30 
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: NEL   ! Nth element
!!
!! Function return value.
      REAL :: EFFECTIVE_STRESS
!!
!! Local variables.
      LOGICAL,        SAVE :: FIRST = .TRUE.
      INTEGER(ITYPE), SAVE :: NHX,NPX,NPY,NTX,NM3,NP3,NM4,NP4,NTR,N
      INTEGER(ITYPE)       :: NDEX(4)
      REAL(RTYPE)          :: TEMP(4)
      REAL(RTYPE)          :: PRESSURE,S1,S2,S3,S4,S5,S6
!!
!! Function reference.
      REAL(RTYPE) :: EFF_STR
!!
      IF (FIRST) THEN
        NHX = NUMHX
        NPX = NHX + NUMPX
        NPY = NPX + NUMPY
        NTX = NPY + NUMTX
        NM3 = NTX + NUMM3
        NP3 = NM3 + NUMP3
        NM4 = NP3 + NUMM4
        NP4 = NM4 + NUMP4
        NTR = NP4 + NUMTR
        FIRST = .FALSE.
      ENDIF
!!
      EFFECTIVE_STRESS = 0.0
      IF (NEL .LE. NHX) THEN
        N = NEL
        PRESSURE = P3RD * SUM( HEXAH(N)%RES%Stress(1:3) )
        S1 = HEXAH(N)%RES%Stress(1) - PRESSURE
        S2 = HEXAH(N)%RES%Stress(2) - PRESSURE
        S3 = HEXAH(N)%RES%Stress(3) - PRESSURE
        S4 = HEXAH(N)%RES%Stress(4) * HEXAH(N)%RES%Stress(4)
        S5 = HEXAH(N)%RES%Stress(5) * HEXAH(N)%RES%Stress(5) 
        S6 = HEXAH(N)%RES%Stress(6) * HEXAH(N)%RES%Stress(6)
        S4 = S4 + S5 + S6
        EFFECTIVE_STRESS = EFF_STR (S1,S2,S3,S4)
      ELSE IF (NEL .LE. NPX) THEN
        N = NEL - NHX
        PRESSURE = P3RD * SUM( PENTA(N)%RES%Stress(1:3) )
        S1 = PENTA(N)%RES%Stress(1) - PRESSURE
        S2 = PENTA(N)%RES%Stress(2) - PRESSURE
        S3 = PENTA(N)%RES%Stress(3) - PRESSURE
        S4 = PENTA(N)%RES%Stress(4) * PENTA(N)%RES%Stress(4)
        S5 = PENTA(N)%RES%Stress(5) * PENTA(N)%RES%Stress(5) 
        S6 = PENTA(N)%RES%Stress(6) * PENTA(N)%RES%Stress(6)
        S4 = S4 + S5 + S6
        EFFECTIVE_STRESS = EFF_STR (S1,S2,S3,S4)
      ELSE IF (NEL .LE. NPY) THEN
        N = NEL - NPX
        PRESSURE = P3RD * SUM( PYRAMID(N)%RES%Stress(1:3) )
        S1 = PYRAMID(N)%RES%Stress(1) - PRESSURE
        S2 = PYRAMID(N)%RES%Stress(2) - PRESSURE
        S3 = PYRAMID(N)%RES%Stress(3) - PRESSURE
        S4 = PYRAMID(N)%RES%Stress(4) * PYRAMID(N)%RES%Stress(4)
        S5 = PYRAMID(N)%RES%Stress(5) * PYRAMID(N)%RES%Stress(5) 
        S6 = PYRAMID(N)%RES%Stress(6) * PYRAMID(N)%RES%Stress(6)
        S4 = S4 + S5 + S6
        EFFECTIVE_STRESS = EFF_STR (S1,S2,S3,S4)
      ELSE IF (NEL .LE. NTX) THEN
        N = NEL - NPY
        Iform = SECTION_3D( TETRA(N)%PAR%SecID )%Iform
        IF (Iform .EQ. Nodal_Based) THEN
          NDEX = TETRA(N)%RES%IN
          TEMP = (/ZERO,ZERO,ZERO,ZERO/)
          TEMP = TEMP + TETNP(NDEX)%Stress(1)
          TEMP = TEMP + TETNP(NDEX)%Stress(2)
          TEMP = TEMP + TETNP(NDEX)%Stress(3)
          PRESSURE = P12TH * SUM( TEMP ) 
        ELSE
          PRESSURE = P3RD * SUM( TETRA(N)%RES%Stress(1:3) )
        ENDIF
        IF (Iform .EQ. Nodal_Based) THEN
          S1 =  P4TH * SUM( TETNP(NDEX)%Stress(1) ) - PRESSURE
          S2 =  P4TH * SUM( TETNP(NDEX)%Stress(2) ) - PRESSURE
          S3 =  P4TH * SUM( TETNP(NDEX)%Stress(3) ) - PRESSURE
          S4 = (P4TH * SUM( TETNP(NDEX)%Stress(4) ) )**2
          S5 = (P4TH * SUM( TETNP(NDEX)%Stress(5) ) )**2
          S6 = (P4TH * SUM( TETNP(NDEX)%Stress(6) ) )**2
          S4 = (S4 + S5 + S6)
        ELSE
          S1 = TETRA(N)%RES%Stress(1) - PRESSURE
          S2 = TETRA(N)%RES%Stress(2) - PRESSURE
          S3 = TETRA(N)%RES%Stress(3) - PRESSURE
          S4 = TETRA(N)%RES%Stress(4) * TETRA(N)%RES%Stress(4)
          S5 = TETRA(N)%RES%Stress(5) * TETRA(N)%RES%Stress(5) 
          S6 = TETRA(N)%RES%Stress(6) * TETRA(N)%RES%Stress(6)
          S4 = S4 + S5 + S6
        ENDIF
        EFFECTIVE_STRESS = EFF_STR (S1,S2,S3,S4)
      ELSE IF (NEL .LE. NM3) THEN
        N = NEL - NTX
        PRESSURE = P3RD * SUM ( MEMBT(N)%RES%Stress(1:2) )
        S1 = MEMBT(N)%RES%Stress(1) - PRESSURE
        S2 = MEMBT(N)%RES%Stress(2) - PRESSURE
        S3 =                        - PRESSURE
        S4 = MEMBT(N)%RES%Stress(3) * MEMBT(N)%RES%Stress(3)
        EFFECTIVE_STRESS = EFF_STR (S1,S2,S3,S4)
      ELSE IF (NEL .LE. NP3) THEN
        N = NEL - NM3
        Ist = PLATT(N)%PAR%Ist
        PRESSURE = P3RD * SUM( STRESS(1:2,Ist) )
        S1 = STRESS(1,Ist) - PRESSURE
        S2 = STRESS(2,Ist) - PRESSURE
        S3 =               - PRESSURE
        S4 = STRESS(4,Ist) * STRESS(4,Ist) 
        EFFECTIVE_STRESS = EFF_STR (S1,S2,S3,S4)
      ELSE IF (NEL .LE. NM4) THEN
        N = NEL - NP3
        PRESSURE = P3RD * SUM( MEMBQ(N)%RES%Stress(1:2) )
        S1 = MEMBQ(N)%RES%Stress(1) - PRESSURE
        S2 = MEMBQ(N)%RES%Stress(2) - PRESSURE
        S3 =                        - PRESSURE
        S4 = MEMBQ(N)%RES%Stress(3) * MEMBQ(N)%RES%Stress(3)
        EFFECTIVE_STRESS = EFF_STR (S1,S2,S3,S4)
      ELSE IF (NEL .LE. NP4) THEN
        N = NEL - NM4 
        Ist = PLATQ(N)%PAR%Ist
        PRESSURE = P3RD * SUM( STRESS(1:2,Ist) )
        S1 = STRESS(1,Ist) - PRESSURE
        S2 = STRESS(2,Ist) - PRESSURE
        S3 =               - PRESSURE
        S4 = STRESS(4,Ist) * STRESS(4,Ist) 
        EFFECTIVE_STRESS = EFF_STR (S1,S2,S3,S4)
      ELSE IF (NEL .LE. NTR) THEN
        N = NEL - NP4
        PRESSURE = P3RD * TRUSS(N)%RES%Stress
        S1 = TRUSS(N)%RES%Stress - PRESSURE
        S2 =                     - PRESSURE
        S3 =                     - PRESSURE
        S4 = 0.0
        EFFECTIVE_STRESS = EFF_STR (S1,S2,S3,S4)
      ENDIF

      RETURN
      END FUNCTION EFFECTIVE_STRESS
!!
      FUNCTION SEGMENTS_AND_NODES_USED ( M )
!!
!! Copyright (c) by FMA Development, LLC, 26-MAY-2004
!!
!! Purpose: Mark segments and nodes used, count number of unique segments
!! and nodes used by the segment set, and create index translation arrays 
!! for segment set "M".
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: M   !  Segment set ID
!!
!! Function return value.
      INTEGER(ITYPE) :: SEGMENTS_AND_NODES_USED
!!
!! Local Variables.
      INTEGER(ITYPE), SAVE :: IX(4),N,K
!!
!! Clear marker/index-sequence and translation arrays.
!!
      NELUSED = 0
      NPTUSED = 0
      NPTNOWI = 0
!!
!! Mark nodes used by this segment set.
!!
      SEGMENTS_AND_NODES_USED = 0

      N = 0
      DO WHILE (NEXT_SEG_ID(M,N))
        NELUSED(N) = 1
        IX(1:4) = SEGMENT(N)%PAR%IX(1:4)
        IF (IX(4) .LE. 0) THEN
          NPTUSED(IX(1:3)) = (/1,1,1/)            
        ELSE
          NPTUSED(IX(1:4)) = (/1,1,1,1/)            
        ENDIF
        SEGMENTS_AND_NODES_USED = 1
      ENDDO
!!
!! Create nodal point index map NPTNOWI for use with 
!! element connectivity n-tuples. (Maps program index 
!! to index written to esg files for this segment set.)
!!
!! Also, covert NPTUSED into a sequential index map.
!!
      K = 0
      DO N = 1,NUMNP
        IF (NPTUSED(N) .EQ. 1) THEN
          K = K + 1
          NPTNOWI(N) = K
          NPTUSED(K) = N
        ENDIF
      ENDDO
      NNodes = K
!!
!! Covert NELUSED into a sequential index map.
!!
      K= 0 
      DO N = 1,NUMSG
        IF (NELUSED(N) .EQ. 1) THEN
          K = K + 1
          NELUSED(K) = N
        ENDIF
      ENDDO
      NElems = K

      RETURN
      END FUNCTION SEGMENTS_AND_NODES_USED
!!
      FUNCTION SEGMENT_N_TUPLE( NEL )
!!
!! Copyright (c) by FMA Development, LLC, 26-MAY-2004
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: NEL   ! Segment ID
!!
!! Function return value.
      INTEGER(ITYPE) :: SEGMENT_N_TUPLE
!!
!! Local variables.
      INTEGER(ITYPE), SAVE :: IX(4)
!!
      IX(1:4) = SEGMENT(NEL)%PAR%IX(1:4)
      IF (IX(4) .LE. 0) THEN
        IX(4) = IX(3)
      ENDIF
      SEGMENT_N_TUPLE = 4
      NTUPLE(1:4) = NPTNOWI(IX(1:4))

      RETURN
      END FUNCTION SEGMENT_N_TUPLE
!!
      FUNCTION WEDGES_AND_NODES_USED ( M )
!!
!! Copyright (c) by FMA Development, LLC, 26-MAY-2004
!!
!! Purpose: Mark wedges and nodes used, count number of unique wedges
!! and nodes used by the wedge set, and create index translation arrays 
!! for wedge set "M".
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: M   ! Solid-to-solid interface ID
!!
!! Function return value.
      INTEGER(ITYPE) :: WEDGES_AND_NODES_USED
!!
!! Local Variables.
      INTEGER(ITYPE), SAVE :: IX(4),N,K
      LOGICAL              :: FOUND
!!
!! Clear marker/index-sequence and translation arrays.
!!
      NELUSED = 0
      NPTUSED = 0
      NPTNOWI = 0
!!
!! Mark wedges used by this solid-to-solid interface.
!!
      WEDGES_AND_NODES_USED = 0

      DO N = 1,NUMWX
        IF (WEDGE(N)%MPCID .EQ. M) THEN
          NELUSED(N) = 1
          WEDGES_AND_NODES_USED = 1
        ENDIF
      ENDDO
!!
!! Note: The wedges don't use the nodal points of the
!! mesh, but contain information for building their
!! own vertex coordinates, see WEDGE_NP_COORDINATE.
!!
!! Create nodal point index map NPTNOWI for use with 
!! element connectivity n-tuples. (Maps program index 
!! to index written to esg files for this wedge set.)
!!
!! Also, covert NPTUSED into a sequential index map.
!!
!!      K = 0
!!      DO N = 1,NUMNP
!!        IF (NPTUSED(N) .EQ. 1) THEN
!!          K = K + 1
!!          NPTNOWI(N) = K
!!          NPTUSED(K) = N
!!        ENDIF
!!      ENDDO
!!      NNodes = K
!!
!! Covert NELUSED into a sequential index map.
!!
      K= 0 
      DO N = 1,NUMWX
        IF (NELUSED(N) .EQ. 1) THEN
          K = K + 1
          NELUSED(K) = N
        ENDIF
      ENDDO

      NElems = K
      NNodes = 6 * NElems

      RETURN
      END FUNCTION WEDGES_AND_NODES_USED
!!
      FUNCTION WEDGE_N_TUPLE( NEL )
!!
!! Copyright (c) by FMA Development, LLC, 26-MAY-2004
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: NEL   ! Wedge ID
!!
!! Function return value.
      INTEGER(ITYPE) :: WEDGE_N_TUPLE
!!
!! This routine supplies element connectivity n-tuples.
!!
      IOS = 6*NEL - 6

      NTUPLE(1:6) = (/1+IOS,2+IOS,3+IOS,4+IOS,5+IOS,6+IOS/)

      WEDGE_N_TUPLE = 6

      RETURN
      END FUNCTION WEDGE_N_TUPLE
!!
      FUNCTION WEDGE_NP_COORDINATE( I, J )
!!
!! Copyright (c) by FMA Development, LLC, 26-MAY-2004
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: I     ! Nodal counter
      INTEGER(ITYPE), INTENT(IN) :: J     ! 1/2/3=x/y/z
!!
!! Function return value.
      REAL(RTYPE) :: WEDGE_NP_COORDINATE
!!
!! Local variables.
      REAL(RTYPE)    :: PXYZ,QXYZ,W
      INTEGER(ITYPE) :: JXYZ,NEL,NPT,N
!!
!! Caution: This routine contains a hard-coded dependency.
!! This routine will generate nodal point coordinates for
!! the sequence of wedges contained in NELUSED(*). Each
!! wedge has its own six nodal points.
!!
      NEL = ((I-1)/6) + 1
      NPT = I - 6*((I-1)/6)

      MEL = NELUSED(NEL)

      PXYZ = ZERO

      DO k = 1,WEDGE(MEL)%K(NPT)
        N =    WEDGE(MEL)%N(k,NPT)
        W =    WEDGE(MEL)%W(k,NPT)
        JXYZ = J
        SELECT CASE (JXYZ)
          CASE(1);  QXYZ = W * MOTION(N)%Px
          CASE(2);  QXYZ = W * MOTION(N)%Py
          CASE(3);  QXYZ = W * MOTION(N)%Pz
        END SELECT
        PXYZ = PXYZ + QXYZ
      ENDDO

      WEDGE_NP_COORDINATE = PXYZ

      RETURN
      END FUNCTION WEDGE_NP_COORDINATE
!!
      FUNCTION WEDGE_NP_VELOCITY( I, J )
!!
!! Copyright (c) by FMA Development, LLC, 26-MAY-2004
!!
!! Arguments.
      INTEGER(ITYPE), INTENT(IN) :: I     ! Nodal counter
      INTEGER(ITYPE), INTENT(IN) :: J     ! 1/2/3=x/y/z
!!
!! Function return value.
      REAL(RTYPE) :: WEDGE_NP_VELOCITY
!!
!! Local variables.
      REAL(RTYPE)    :: VXYZ,QXYZ,W
      INTEGER(ITYPE) :: JXYZ,NEL,NPT,N
!!
!! Caution: This routine contains a hard-coded dependency.
!! This routine will generate nodal point velocities for
!! the sequence of wedges contained in NELUSED(*). Each
!! wedge has its own six nodal points.
!!
      NEL = ((I-1)/6) + 1
      NPT = I - 6*((I-1)/6)

      MEL = NELUSED(NEL)

      VXYZ = ZERO

      DO k = 1,WEDGE(MEL)%K(NPT)
        N =    WEDGE(MEL)%N(k,NPT)
        W =    WEDGE(MEL)%W(k,NPT)
        JXYZ = J
        SELECT CASE (JXYZ)
          CASE(1);  QXYZ = W * MOTION(N)%Vx
          CASE(2);  QXYZ = W * MOTION(N)%Vy
          CASE(3);  QXYZ = W * MOTION(N)%Vz
        END SELECT
        VXYZ = VXYZ + QXYZ
      ENDDO

      WEDGE_NP_VELOCITY = VXYZ

      RETURN
      END FUNCTION WEDGE_NP_VELOCITY

      END SUBROUTINE WRITE_TO_ESG_RESULTS_FILE


More information about the ParaView mailing list