c$Id: umacr9.f,v 1.1 2008/12/12 23:58:47 rlt Exp $ subroutine umacr9(lct,ctl) c * * F E A P * * A Finite Element Analysis Program c.... Copyright (c) 1984-2009: Regents of the University of California c All rights reserved c-----[--.----+----.----+----.-----------------------------------------] c Modification log Date (dd/mm/year) c Original version 01/11/2006 c 1. Remove 'prt' from argument list 09/07/2009 c-----[--.----+----.----+----.-----------------------------------------] c Purpose: User interface for adding solution command language c instructions. c Inputs: c lct - Command character parameters c ctl(3) - Command numerical parameters c Outputs: c N.B. Users are responsible for command actions. See c programmers manual for example. c-----[--.----+----.----+----.-----------------------------------------] implicit none include 'cdata.h' include 'comfil.h' include 'counts.h' include 'iofile.h' include 'pointer.h' include 'sdata.h' include 'strnum.h' include 'umac1.h' include 'comblk.h' logical pcomp character lct*15 character fname*128 integer i,ii,node real*8 ctl(3) save c Set command word if(pcomp(uct,'mac9',4)) then ! Usual form uct = 'pvie' ! Specify 'name' elseif(urest.eq.1) then ! Read restart data elseif(urest.eq.2) then ! Write restart data else ! Perform user operation if(.not.pcomp(lct,'time',4)) then open(unit=99,file=lct,access='sequential') write(*,*) 'Saving PARAVIEW data to ',lct else i = index(fplt,' ') fname(1:128) = ' ' fname(1:i-1)=fplt(1:i-1) fname(i:i+4) = '00000' if (nstep.le.9) then write(fname(i+4:i+4),'(i1)') nstep elseif (nstep.le.99) then write(fname(i+3:i+4),'(i2)') nstep elseif (nstep.le.999) then write(fname(i+2:i+4),'(i3)') nstep elseif (nstep.le.9999) then write(fname(i+1:i+4),'(i4)') nstep elseif (nstep.le.99999) then write(fname(i:i+4),'(i5)') nstep endif call addext(fname,'vtu',128,3) open(unit=99,file=fname,access='sequential') write(*,*) 'Saving PARAVIEW data to ',fname end if c Write out top header write(99,1000) write(99,1010) numnp,numel ! Start Mesh/Piece Section write(99,1020) ! Start Point/Node data write(99,1030) 'Float64','nodes',ndm do i = 1,numnp do ii = 1,ndm write(99,5000) hr(np(43)+(i-1)*ndm+(ii-1)) end do end do write(99,1031) ! Close Node data write(99,1021) ! Close Points section write(99,1050) ! Start Cell Section write(99,1030) 'Int32','connectivity',1 ! Start Elements c Offsets memory allocation call palloc(111,'TEMP1',numel+1,1) mr(np(111)) = 0; do i = 1,numel mr(np(111)+i) = mr(np(111)+i-1) do ii = 1,nen node = mr(np(33) + ii-1 + nen1*(i-1)) if (node .ne. 0) then write(99,6000) node-1 mr(np(111)+i) = mr(np(111)+i) + 1 endif end do end do write(99,1031) ! Close Elements write(99,1030) 'Int32','offsets',1 ! Start Offsets do i = 1,numel write(99,6000) mr(np(111)+i) end do write(99,1031) ! Close Offsets write(99,1030) 'UInt8','types',1 ! Start Element types do i = 1,numel if (mr(np(111)+i)-mr(np(111)+i-1) .eq. 2) then ! 2 node line write(99,6000) 3 elseif (mr(np(111)+i)-mr(np(111)+i-1) .eq. 3) then ! 3 node triangle write(99,6000) 5 elseif (mr(np(111)+i)-mr(np(111)+i-1) .eq. 4) then ! 4 node quad write(99,6000) 9 elseif (mr(np(111)+i)-mr(np(111)+i-1) .eq. 8) then ! 8 node brick write(99,6000) 12 elseif (mr(np(111)+i)-mr(np(111)+i-1) .eq. 27) then ! 27 node brick (plot as 20 node) write(99,6000) 25 endif end do call palloc(111,'TEMP1',0,1) write(99,1031) ! Close Element types write(99,1051) ! Close Cell Section write(99,1060) ! Start Point Data write(99,1030) 'Float64','dist',ndf ! Start Displacements do i = 1,numnp do ii =0,ndf-1 write(99,5000) hr(np(40)+(i-1)*ndf+ii) end do end do write(99,1031) ! Close Displacements if (istv.ne.0) then write(99,1030) 'Float64','Stresses',abs(istv) ! Start Stresses do i = 1,numnp do ii =1,abs(istv) write(99,5000) hr(np(58)+(i-1) + ii*numnp) end do end do write(99,1031) ! Close Stresses write(99,1030) 'Float64','PStress',7 ! Start Principal Stresses do i = 1,numnp do ii =1,7 write(99,5000) hr(np(57)+(i-1) + ii*numnp) end do end do write(99,1031) ! Close Stresses endif write(99,1061) ! Close Point Data Section write(99,1011) ! Close Mesh/Piece c Close the XML file write(99,1001) close(99) endif 1000 format('',/ * '',/ * '') 1001 format(' ') 1010 format('') 1011 format('') 1020 format('') 1021 format('') 1030 format('') 1031 format('') 1060 format('') 1061 format('') c 1070 format('') c 1071 format('') 1050 format('') 1051 format('') 5000 format(e15.5,' ',$) 6000 format(i6,' ',$) end