program f90_FEM_FRAME !******************************************************* !2D frame analysis !******************************************************* implicit none integer::i,j,k,l,m,kk,ne,n,ia,ja character::strcom*50 !Comment integer::NODT !Number of nodes integer::NELT !Number of elements integer::MATEL !Number of material sets integer::KOX !Number of restricted nodes in x-direction integer::KOY !Number of restricted nodes in y-direction integer::KOZ !Number of restricted nodes in rotation integer::NF !Number of loaded nodes integer,allocatable::kakom(:,:) !Element connectivity integer,allocatable::matno(:) !Material set number real(8),allocatable::wm(:,:) !work for material properties real(8),allocatable::Em(:) !Elastic modulus of element real(8),allocatable::AA(:) !Section area of element real(8),allocatable::AI(:) !Moment of inertia real(8),allocatable::gamma(:) !Unit weight of element real(8),allocatable::gkh(:) !Horizontal acceleration for element in the ratio to 'g' real(8),allocatable::gkv(:) !Vertical acceleration for element in the ratio to 'g' real(8),allocatable::alpha(:) !Coefficient of thermal expansion of element real(8),allocatable::x(:) !x-coordinate of node real(8),allocatable::y(:) !y-coordinate of node real(8),allocatable::deltaT(:) !Temperature change of node (+: temperature zise) integer,allocatable::nokx(:) !Restricted node number in x-direction integer,allocatable::noky(:) !Restricted node number in y-direction integer,allocatable::nokz(:) !Restricted node number in rotation direction real(8),allocatable::f(:) !External force vector integer,allocatable::mindex(:) !Index for treatment of restricted nodes real(8),allocatable::ftvec(:) !Load vector real(8),allocatable::ftemp(:) !Load vector due to temperature change real(8),allocatable::fmass(:) !Load vector due to body force real(8),allocatable::fdis(:) !Load vector caused by restricted nodal displacements real(8),allocatable::rdis(:) !Displacement of restricted nodes real(8),allocatable::tsm(:,:) !Stiffness matrix real(8),allocatable::dis(:) !Nodal displacement vector real(8),allocatable::reac(:) !work vector real(8),allocatable::sm(:,:) !Element stiffness matrix real(8),allocatable::hen(:,:) !Coordinate transform matrix real(8),allocatable::fevec(:) !Element nodal force matrix integer::nt !Number of degrees of freedom of all nodes (NODT x nhen) integer::nod=2 !Number of nodes per element integer::nhen=3 !Number of degrees of freedom per node integer::mm !Number of degrees of freedom of reduced FE equation (NODT x nhen) integer::ib !Band width of reduced matrix real(8)::fact0=1.0D-30 !zero judgement character::fnameR*50,fnameW*50 !Input file name, output file name character :: linebuf*1000 character :: fmt1*200,fmt2*200,fmt3*200,fmt4*200 real(4) :: t1,t2 fmt1="(i5,2(',',e15.7),3(',',e15.7),3(',',i2),3(',',e15.7),',',e15.7)" !Nodes fmt2="(i5,2(',',i5),7(',',e15.7),',',i5)" !Elements fmt3="(i5,2(',',e15.7),3(',',e15.7),3(',',e15.7),3(',',e15.7))" !Displacement & reaction fmt4="(i5,3(',',e15.7),3(',',e15.7))" !Stress resultant !To input file names call getarg(1,fnameR) !Input file name call getarg(2,fnameW) !Output file name !To start the measurement of calculation time t1=DIFT() !*************************** !data input !*************************** open(11,file=fnameR,status='old') !/*--------------------------------------------------------------------------------------*/ !/*Comment*/ !/*--------------------------------------------------------------------------------------*/ read(11,'(a)') strcom !------------------------------------------------------------------------------------------ !Basic values for analysis !------------------------------------------------------------------------------------------ read(11,*) NODT,NELT,MATEL,KOX,KOY,KOZ,NF !------------------------------------------------------------------------------------------ !To declare array size !------------------------------------------------------------------------------------------ nt=NODT*nhen allocate(kakom(1:NELT,1:nod)) allocate(matno(1:NELT)) allocate(wm(1:NELT,1:7)) allocate(Em(1:NELT)) allocate(AA(1:NELT)) allocate(AI(1:NELT)) allocate(gamma(1:NELT)) allocate(gkh(1:NELT)) allocate(gkv(1:NELT)) allocate(alpha(1:NELT)) allocate(x(1:NODT)) allocate(y(1:NODT)) allocate(deltaT(1:NODT)) if(1<=KOX)allocate(nokx(1:KOX)) if(1<=KOY)allocate(noky(1:KOY)) if(1<=KOZ)allocate(nokz(1:KOZ)) allocate(f(1:nt)) allocate(mindex(1:nt)) allocate(ftvec(1:nt)) allocate(ftemp(1:nt)) allocate(fmass(1:nt)) allocate(fdis(1:nt)) allocate(rdis(1:nt)) allocate(tsm(1:nt,1:nt)) allocate(dis(1:nt)) allocate(reac(1:nt)) allocate(sm(1:nod*nhen,1:nod*nhen)) allocate(hen(1:nod*nhen,1:nod*nhen)) allocate(fevec(1:nod*nhen)) !------------------------------------------------------------------------------------------ !Material properties (Em,AA,AI,gamma,gkh,gkv,alpha) !------------------------------------------------------------------------------------------ do i=1,MATEL read(11,*) wm(i,1),wm(i,2),wm(i,3),wm(i,4),wm(i,5),wm(i,6),wm(i,7) end do !------------------------------------------------------------------------------------------ !Element-nodes relationship, material set No. !------------------------------------------------------------------------------------------ do ne=1,NELT read(11,*) kakom(ne,1),kakom(ne,2),matno(ne) do i=1,MATEL if(i==matno(ne))then Em(ne) =wm(i,1) AA(ne) =wm(i,2) AI(ne) =wm(i,3) gamma(ne)=wm(i,4) gkh(ne) =wm(i,5) gkv(ne) =wm(i,6) alpha(ne)=wm(i,7) end if end do end do deallocate(wm) !------------------------------------------------------------------------------------------ !(x,y) coordinates of node, temperature change of node !------------------------------------------------------------------------------------------ do i=1,NODT read(11,*) x(i),y(i),deltaT(i) end do !------------------------------------------------------------------------------------------ !Restricted node number, displacement of restricted node !------------------------------------------------------------------------------------------ do i=1,nt rdis(i)=0.0D0 end do if(00)then mm=mm+1 mindex(mm)=i end if end do do i=1,mm ia=mindex(i) reac(i)=reac(ia) do j=1,mm ja=mindex(j) tsm(i,j)=tsm(ia,ja) end do end do !----------------------------------------------------- !Calculation of band width and memory of band matrix !----------------------------------------------------- ib=IBAND(mm,tsm,fact0,nt) call BAND(mm,ib,tsm,nt) !----------------------------------------------------- !To verify diagonal element & to stop if necessary !----------------------------------------------------- do i=1,mm if(abs(tsm(i,1))m1)then mm=m1 else mm=n1-i+1 end if do j=1,mm z=array(i,j) if(j/=m1)then if(m1-j+11)then array(i,j)=z/array(i,1) else array(i,1)=sqrt(z) end if else array(i,j)=z/array(i,1) end if end do end do end subroutine CHBAND10 subroutine CHBAND11(n1,m1,array,vector,nt) integer,intent(in)::n1,m1,nt real(8),intent(in)::array(1:nt,1:nt) real(8),intent(out)::vector(nt) integer :: i,j,mm real(8) :: z vector(1)=vector(1)/array(1,1) do i=2,n1 z=vector(i) if(i>m1)then mm=m1 else mm=i end if do j=2,mm z=z-array(i-j+1,j)*vector(i-j+1) end do vector(i)=z/array(i,1) end do vector(n1)=vector(n1)/array(n1,1) do i=n1-1,1,-1 z=vector(i) do j=2,m1 if(i+j-1>n1)exit z=z-array(i,j)*vector(i+j-1) end do vector(i)=z/array(i,1) end do end subroutine CHBAND11 integer function IBAND(mm,array,fact0,nt) integer,intent(in)::mm,nt real(8),intent(in)::array(1:nt,1:nt),fact0 integer :: i,j !----------------------------------------------------- !calculation of band width !----------------------------------------------------- IBAND=1 do i=1,mm-1 do j=mm,i+1,-1 if(fact0<=abs(array(i,j)))exit end do if(IBAND<=j-i+1)IBAND=j-i+1 end do if(mm