! ! CalculiX - A 3-dimensional finite element program ! Copyright (C) 1998-2007 Guido Dhondt ! ! This program is free software; you can redistribute it and/or ! modify it under the terms of the GNU General Public License as ! published by the Free Software Foundation(version 2); ! ! ! This program is distributed in the hope that it will be useful, ! but WITHOUT ANY WARRANTY; without even the implied warranty of ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the ! GNU General Public License for more details. ! ! You should have received a copy of the GNU General Public License ! along with this program; if not, write to the Free Software ! Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA. ! subroutine mafillsmcs(co,nk,kon,ipkon,lakon,ne,nodeboun,ndirboun, & xboun,nboun, & ipompc,nodempc,coefmpc,nmpc,nodeforc,ndirforc,xforc, & nforc,nelemload,sideload,xload,nload,xbody,ipobody,nbody,cgr, & ad,au,bb,nactdof,icol,jq,irow,neq,nzl,nmethod, & ikmpc,ilmpc,ikboun,ilboun,elcon,nelcon,rhcon, & nrhcon,alcon,nalcon,alzero,ielmat,ielorien,norien,orab,ntmat_, & t0,t1,ithermal,prestr, & iprestr,vold,iperturb,sti,nzs,stx,adb,aub,iexpl,plicon, & nplicon,plkcon,nplkcon,xstiff,npmat_,dtime, & matname,mint_,ics,cs,nm,ncmat_,labmpc,mass,stiffness,buckling, & rhsi,intscheme,mcs,coriolis,ibody) ! ! filling the stiffness matrix in spare matrix format (sm) ! for cyclic symmetry calculations ! implicit none ! logical mass,stiffness,buckling,rhsi,coriolis ! character*8 lakon(*) character*20 labmpc(*),sideload(*) character*80 matname(*) ! integer kon(*),nodeboun(*),ndirboun(*),ipompc(*),nodempc(3,*), & nodeforc(2,*),ndirforc(*),nelemload(2,*),icol(*),jq(*),ikmpc(*), & ilmpc(*),ikboun(*),ilboun(*),nactdof(0:3,*),konl(20),irow(*), & nelcon(2,*),nrhcon(*),nalcon(2,*),ielmat(*),ielorien(*), & ipkon(*),ics(*),ij,ilength,lprev,ipobody(2,*),nbody, & ibody(3,*) ! integer nk,ne,nboun,nmpc,nforc,nload,neq,nzl,nmethod, & ithermal,iprestr,iperturb,nzs,i,j,k,l,m,idist,jj, & ll,id,id1,id2,ist,ist1,ist2,index,jdof1,jdof2,idof1,idof2, & mpc1,mpc2,index1,index2,node1,node2,kflag, & ntmat_,indexe,nope,norien,iexpl,mint_,i0,nm,inode,icomplex, & inode1,icomplex1,inode2,icomplex2,ner,ncmat_,intscheme,istep, & iinc,mcs ! integer nplicon(0:ntmat_,*),nplkcon(0:ntmat_,*),npmat_ ! real*8 co(3,*),xboun(*),coefmpc(*),xforc(*),xload(2,*),p1(3), & p2(3),ad(*),au(*),bodyf(3),bb(*),xbody(7,*),cgr(4,*), & t0(*),t1(*),prestr(6,mint_,*),vold(0:3,*),s(60,60),ff(60), & sti(6,mint_,*),sm(60,60),stx(6,mint_,*),adb(*),aub(*), & elcon(0:ncmat_,ntmat_,*),rhcon(0:1,ntmat_,*), & alcon(0:6,ntmat_,*),cs(17,*),alzero(*),orab(7,*) ! real*8 plicon(0:2*npmat_,ntmat_,*),plkcon(0:2*npmat_,ntmat_,*), & xstiff(27,mint_,*),pi,theta,ti,tr ! real*8 om,valu2,value,dtime,walue,walu2,time,ttime ! ! ! calculating the scaling factors for the cyclic symmetry calculation ! c do i=1,nmpc c write(*,*) i,labmpc(i) c index=ipompc(i) c do c write(*,'(i5,1x,i5,1x,e11.4)') nodempc(1,index), c & nodempc(2,index),coefmpc(index) c index=nodempc(3,index) c if(index.eq.0) exit c enddo c write(*,*) c enddo ! pi=4.d0*datan(1.d0) ! do i=1,mcs theta=nm*2.d0*pi/cs(1,i) cs(15,i)=dcos(theta) cs(16,i)=dsin(theta) enddo ! kflag=2 i0=0 ! ! determining nzl ! nzl=0 do i=neq,1,-1 if(icol(i).gt.0) then nzl=i exit endif enddo ! ! initializing the matrices ! do i=1,neq ad(i)=0.d0 enddo do i=1,nzs au(i)=0.d0 enddo ! do i=1,neq adb(i)=0.d0 enddo do i=1,nzs aub(i)=0.d0 enddo ! ner=neq/2 ! ! loop over all elements ! ! initialisation of the error parameter ! do i=1,ne ! if(ipkon(i).lt.0) cycle indexe=ipkon(i) if(lakon(i)(4:4).eq.'2') then nope=20 elseif(lakon(i)(4:4).eq.'8') then nope=8 elseif(lakon(i)(4:5).eq.'10') then nope=10 elseif(lakon(i)(4:4).eq.'4') then nope=4 elseif(lakon(i)(4:5).eq.'15') then nope=15 else nope=6 endif ! do j=1,nope konl(j)=kon(indexe+j) enddo c! c! assigning centrifugal forces c! c j=ipobody(1,i) c if(j.ne.0) then c om=xbody(1,j) c p1(1)=xbody(2,j) c p1(2)=xbody(3,j) c p1(3)=xbody(4,j) c p2(1)=xbody(5,j) c p2(2)=xbody(6,j) c p2(3)=xbody(7,j) c else c om=0.d0 c endif ! om=0.d0 ! if(nbody.gt.0) then ! ! assigning centrifugal forces ! index=i do j=ipobody(1,index) if(j.eq.0) exit if(ibody(1,j).eq.1) then om=xbody(1,j) p1(1)=xbody(2,j) p1(2)=xbody(3,j) p1(3)=xbody(4,j) p2(1)=xbody(5,j) p2(2)=xbody(6,j) p2(3)=xbody(7,j) endif index=ipobody(2,index) if(index.eq.0) exit enddo endif ! call e_c3d(co,nk,konl,lakon(i),p1,p2,om,bodyf,nbody,s,sm,ff,i, & nmethod,elcon,nelcon,rhcon,nrhcon,alcon,nalcon, & alzero,ielmat,ielorien,norien,orab,ntmat_, & t0,t1,ithermal,vold,iperturb, & nelemload,sideload,xload,nload,idist,sti,stx, & iexpl,plicon, & nplicon,plkcon,nplkcon,xstiff,npmat_, & dtime,matname,mint_,ncmat_,mass,stiffness,buckling,rhsi, & intscheme,ttime,time,istep,iinc,coriolis) ! do jj=1,3*nope ! j=(jj-1)/3+1 k=jj-3*(j-1) ! node1=kon(indexe+j) jdof1=nactdof(k,node1) ! do ll=jj,3*nope ! l=(ll-1)/3+1 m=ll-3*(l-1) ! node2=kon(indexe+l) jdof2=nactdof(m,node2) ! ! check whether one of the DOF belongs to a SPC or MPC ! if((jdof1.ne.0).and.(jdof2.ne.0)) then call add_sm_ei(au,ad,aub,adb,jq,irow,jdof1,jdof2, & s(jj,ll),sm(jj,ll),jj,ll) call add_sm_ei(au,ad,aub,adb,jq,irow,jdof1+ner,jdof2+ner, & s(jj,ll),sm(jj,ll),jj,ll) elseif((jdof1.ne.0).or.(jdof2.ne.0)) then ! ! idof1: genuine DOF ! idof2: nominal DOF of the SPC/MPC ! if(jdof1.eq.0) then idof1=jdof2 idof2=(node1-1)*7+k else idof1=jdof1 idof2=(node2-1)*7+m endif ! if(nmpc.gt.0) then call nident(ikmpc,idof2,nmpc,id) if((id.gt.0).and.(ikmpc(id).eq.idof2)) then ! ! regular DOF / MPC ! id1=ilmpc(id) ist=ipompc(id1) index=nodempc(3,ist) if(index.eq.0) cycle do inode=nodempc(1,index) icomplex=0 c write(*,*) id1,labmpc(id1)(1:9) if(labmpc(id1)(1:6).eq.'CYCLIC') then read(labmpc(id1)(7:20),'(i14)') icomplex elseif(labmpc(id1)(1:9).eq.'SUBCYCLIC') then do ij=1,mcs ilength=int(cs(4,ij)) lprev=int(cs(14,ij)) call nident(ics(lprev+1),inode,ilength,id) if(id.gt.0) then if(ics(lprev+id).eq.inode) then icomplex=ij exit endif endif enddo endif idof2=nactdof(nodempc(2,index),inode) if(idof2.ne.0) then value=-coefmpc(index)*s(jj,ll)/coefmpc(ist) valu2=-coefmpc(index)*sm(jj,ll)/ & coefmpc(ist) if(idof1.eq.idof2) then value=2.d0*value valu2=2.d0*valu2 endif if(icomplex.eq.0) then call add_sm_ei(au,ad,aub,adb,jq,irow, & idof1,idof2,value,valu2,i0,i0) call add_sm_ei(au,ad,aub,adb,jq,irow, & idof1+ner,idof2+ner,value,valu2,i0,i0) else walue=value*cs(15,icomplex) walu2=valu2*cs(15,icomplex) call add_sm_ei(au,ad,aub,adb,jq,irow, & idof1,idof2,walue,walu2,i0,i0) call add_sm_ei(au,ad,aub,adb,jq,irow, & idof1+ner,idof2+ner,walue,walu2,i0,i0) if(idof1.ne.idof2) then walue=value*cs(16,icomplex) walu2=valu2*cs(16,icomplex) call add_sm_ei(au,ad,aub,adb,jq,irow, & idof1,idof2+ner,walue,walu2,i0,i0) walue=-walue walu2=-walu2 call add_sm_ei(au,ad,aub,adb,jq,irow, & idof1+ner,idof2,walue,walu2,i0,i0) endif endif endif index=nodempc(3,index) if(index.eq.0) exit enddo cycle endif endif ! else idof1=(node1-1)*7+k idof2=(node2-1)*7+m ! mpc1=0 mpc2=0 if(nmpc.gt.0) then call nident(ikmpc,idof1,nmpc,id1) if((id1.gt.0).and.(ikmpc(id1).eq.idof1)) mpc1=1 call nident(ikmpc,idof2,nmpc,id2) if((id2.gt.0).and.(ikmpc(id2).eq.idof2)) mpc2=1 endif if((mpc1.eq.1).and.(mpc2.eq.1)) then id1=ilmpc(id1) id2=ilmpc(id2) if(id1.eq.id2) then ! ! MPC id1 / MPC id1 ! ist=ipompc(id1) index1=nodempc(3,ist) if(index1.eq.0) cycle do inode1=nodempc(1,index1) icomplex1=0 if(labmpc(id1)(1:6).eq.'CYCLIC') then read(labmpc(id1)(7:20),'(i14)') icomplex1 elseif(labmpc(id1)(1:9).eq.'SUBCYCLIC') then do ij=1,mcs ilength=int(cs(4,ij)) lprev=int(cs(14,ij)) call nident(ics(lprev+1),inode1, & ilength,id) if(id.gt.0) then if(ics(lprev+id).eq.inode1) then icomplex1=ij exit endif endif enddo endif idof1=nactdof(nodempc(2,index1),inode1) index2=index1 do inode2=nodempc(1,index2) icomplex2=0 if(labmpc(id1)(1:6).eq.'CYCLIC') then read(labmpc(id1)(7:20),'(i14)') icomplex2 elseif(labmpc(id1)(1:9).eq.'SUBCYCLIC') then do ij=1,mcs ilength=int(cs(4,ij)) lprev=int(cs(14,ij)) call nident(ics(lprev+1),inode2, & ilength,id) if(id.gt.0) then if(ics(lprev+id).eq.inode2) then icomplex2=ij exit endif endif enddo endif idof2=nactdof(nodempc(2,index2),inode2) if((idof1.ne.0).and.(idof2.ne.0)) then value=coefmpc(index1)*coefmpc(index2)* & s(jj,ll)/coefmpc(ist)/coefmpc(ist) valu2=coefmpc(index1)*coefmpc(index2)* & sm(jj,ll)/coefmpc(ist)/coefmpc(ist) if((icomplex1.eq.0).and. & (icomplex2.eq.0)) then call add_sm_ei(au,ad,aub,adb,jq, & irow,idof1,idof2,value,valu2,i0,i0) call add_sm_ei(au,ad,aub,adb,jq, & irow,idof1+ner,idof2+ner,value, & valu2,i0,i0) elseif((icomplex1.ne.0).and. & (icomplex2.ne.0)) then if(icomplex1.eq.icomplex2) then call add_sm_ei(au,ad,aub,adb,jq, & irow,idof1,idof2,value,valu2,i0,i0) call add_sm_ei(au,ad,aub,adb,jq, & irow,idof1+ner,idof2+ner,value, & valu2,i0,i0) else tr=cs(15,icomplex1)*cs(15,icomplex2) & +cs(16,icomplex1)*cs(16,icomplex2) c write(*,*) 'tr= ',tr walue=value*tr walu2=valu2*tr call add_sm_ei(au,ad,aub,adb,jq, & irow,idof1,idof2,walue,walu2,i0,i0) call add_sm_ei(au,ad,aub,adb,jq, & irow,idof1+ner,idof2+ner,walue, & walu2,i0,i0) ti=cs(15,icomplex1)*cs(16,icomplex2) & -cs(15,icomplex2)*cs(16,icomplex1) c write(*,*) icomplex1,icomplex2, c & cs(15,icomplex1),cs(16,icomplex1), c & cs(15,icomplex2),cs(16,icomplex2) c write(*,*) 'ti= ',ti walue=value*ti walu2=valu2*ti c write(*,'(2i8,2(1x,e11.4))') c & idof1,idof2+ner, c & walue,walu2 call add_sm_ei(au,ad,aub,adb,jq,irow & ,idof1,idof2+ner,walue,walu2,i0,i0) walue=-walue walu2=-walu2 call add_sm_ei(au,ad,aub,adb,jq,irow & ,idof1+ner,idof2,walue,walu2,i0,i0) endif elseif((icomplex1.eq.0).or. & (icomplex2.eq.0)) then if(icomplex2.ne.0) then walue=value*cs(15,icomplex2) walu2=valu2*cs(15,icomplex2) else walue=value*cs(15,icomplex1) walu2=valu2*cs(15,icomplex1) endif call add_sm_ei(au,ad,aub,adb,jq,irow, & idof1,idof2,walue,walu2,i0,i0) call add_sm_ei(au,ad,aub,adb,jq,irow, & idof1+ner,idof2+ner,walue,walu2,i0,i0) if(icomplex2.ne.0) then walue=value*cs(16,icomplex2) walu2=valu2*cs(16,icomplex2) else walue=-value*cs(16,icomplex1) walu2=-valu2*cs(16,icomplex1) endif c walue=value*st c walu2=valu2*st call add_sm_ei(au,ad,aub,adb,jq,irow, & idof1,idof2+ner,walue,walu2,i0,i0) walue=-walue walu2=-walu2 call add_sm_ei(au,ad,aub,adb,jq,irow, & idof1+ner,idof2,walue,walu2,i0,i0) endif endif index2=nodempc(3,index2) if(index2.eq.0) exit enddo index1=nodempc(3,index1) if(index1.eq.0) exit enddo else ! ! MPC id1 / MPC id2 ! ist1=ipompc(id1) index1=nodempc(3,ist1) if(index1.eq.0) cycle do inode1=nodempc(1,index1) icomplex1=0 if(labmpc(id1)(1:6).eq.'CYCLIC') then read(labmpc(id1)(7:20),'(i14)') icomplex1 elseif(labmpc(id1)(1:9).eq.'SUBCYCLIC') then do ij=1,mcs ilength=int(cs(4,ij)) lprev=int(cs(14,ij)) call nident(ics(lprev+1),inode1, & ilength,id) if(id.gt.0) then if(ics(lprev+id).eq.inode1) then icomplex1=ij exit endif endif enddo endif idof1=nactdof(nodempc(2,index1),inode1) ist2=ipompc(id2) index2=nodempc(3,ist2) if(index2.eq.0) then index1=nodempc(3,index1) if(index1.eq.0) then exit else cycle endif endif do inode2=nodempc(1,index2) icomplex2=0 if(labmpc(id2)(1:6).eq.'CYCLIC') then read(labmpc(id2)(7:20),'(i14)') icomplex2 elseif(labmpc(id2)(1:9).eq.'SUBCYCLIC') then do ij=1,mcs ilength=int(cs(4,ij)) lprev=int(cs(14,ij)) call nident(ics(lprev+1),inode2, & ilength,id) if(id.gt.0) then if(ics(lprev+id).eq.inode2) then icomplex2=ij exit endif endif enddo endif idof2=nactdof(nodempc(2,index2),inode2) if((idof1.ne.0).and.(idof2.ne.0)) then value=coefmpc(index1)*coefmpc(index2)* & s(jj,ll)/coefmpc(ist1)/coefmpc(ist2) valu2=coefmpc(index1)*coefmpc(index2)* & sm(jj,ll)/coefmpc(ist1)/coefmpc(ist2) if(idof1.eq.idof2) then value=2.d0*value valu2=2.d0*valu2 endif if((icomplex1.eq.0).and. & (icomplex2.eq.0)) then call add_sm_ei(au,ad,aub,adb,jq, & irow,idof1,idof2,value,valu2,i0,i0) call add_sm_ei(au,ad,aub,adb,jq, & irow,idof1+ner,idof2+ner,value, & valu2,i0,i0) elseif((icomplex1.ne.0).and. & (icomplex2.ne.0)) then if(icomplex1.eq.icomplex2) then call add_sm_ei(au,ad,aub,adb,jq, & irow,idof1,idof2,value,valu2,i0,i0) call add_sm_ei(au,ad,aub,adb,jq, & irow,idof1+ner,idof2+ner,value, & valu2,i0,i0) else tr=cs(15,icomplex1)*cs(15,icomplex2) & +cs(16,icomplex1)*cs(16,icomplex2) c write(*,*) 'tr= ',tr walue=value*tr walu2=valu2*tr call add_sm_ei(au,ad,aub,adb,jq, & irow,idof1,idof2,walue,walu2,i0,i0) call add_sm_ei(au,ad,aub,adb,jq, & irow,idof1+ner,idof2+ner,walue, & walu2,i0,i0) ti=cs(15,icomplex1)*cs(16,icomplex2) & -cs(15,icomplex2)*cs(16,icomplex1) c write(*,*) icomplex1,icomplex2, c & cs(15,icomplex1),cs(16,icomplex1), c & cs(15,icomplex2),cs(16,icomplex2) c write(*,*) 'ti= ',ti walue=value*ti walu2=valu2*ti c write(*,'(2i8,2(1x,e11.4))') c & idof1,idof2+ner, c & walue,walu2 call add_sm_ei(au,ad,aub,adb,jq,irow & ,idof1,idof2+ner,walue,walu2,i0,i0) walue=-walue walu2=-walu2 call add_sm_ei(au,ad,aub,adb,jq,irow & ,idof1+ner,idof2,walue,walu2,i0,i0) endif elseif((icomplex1.eq.0).or. & (icomplex2.eq.0)) then if(icomplex2.ne.0) then walue=value*cs(15,icomplex2) walu2=valu2*cs(15,icomplex2) else walue=value*cs(15,icomplex1) walu2=valu2*cs(15,icomplex1) endif call add_sm_ei(au,ad,aub,adb,jq,irow, & idof1,idof2,walue,walu2,i0,i0) call add_sm_ei(au,ad,aub,adb,jq,irow, & idof1+ner,idof2+ner,walue,walu2,i0,i0) if(idof1.ne.idof2) then if(icomplex2.ne.0) then walue=value*cs(16,icomplex2) walu2=valu2*cs(16,icomplex2) else walue=-value*cs(16,icomplex1) walu2=-valu2*cs(16,icomplex1) endif c walue=value*st c walu2=valu2*st call add_sm_ei(au,ad,aub,adb,jq, & irow,idof1,idof2+ner,walue, & walu2,i0,i0) walue=-walue walu2=-walu2 call add_sm_ei(au,ad,aub,adb,jq, & irow,idof1+ner,idof2,walue, & walu2,i0,i0) endif endif endif index2=nodempc(3,index2) if(index2.eq.0) exit enddo index1=nodempc(3,index1) if(index1.eq.0) exit enddo endif endif endif enddo ! enddo enddo ! return end