%-------------------------- % @Author: Jingqiao Hu % @Date: 2022-09-14 15:04:09 % @LastEditTime: 2022-09-14 15:04:10 %-------------------------- %-------------------------- % @Author: Jingqiao Hu % @Date: 2021-12-09 19:58:02 % @LastEditTime: 2022-07-22 16:03:56 % x = Material indicator matrix. Size used to determine nelx/nely/nelz % dx = the Unit cell length %-------------------------- function Q = homogenization3D_energy_based(nelx, nely, nelz, E0, Emin, nu, penal, x, dx) nelx = 20; nely = 20; nelz = 20; E0 = 1; Emin = 1e-9; nu = 0.3; penal = 3; x = ones(nely, nelx, nelz); dx = 1; nele = nelx*nely*nelz; D0 = E0/(1+nu)/(1-2*nu)* [ 1-nu nu nu 0 0 0; nu 1-nu nu 0 0 0; nu nu 1-nu 0 0 0; 0 0 0 (1-2*nu)/2 0 0; 0 0 0 0 (1-2*nu)/2 0; 0 0 0 0 0 (1-2*nu)/2;]; KE = elementMatVec3D(dx, dx, dx, D0); Num_node = (1+nely)*(1+nelx)*(1+nelz); nodenrs = reshape(1:Num_node,1+nely,1+nelx,1+nelz); edofVec = reshape(3*nodenrs(1:end-1,1:end-1,1:end-1)+1,nelx*nely*nelz,1); edofMat = repmat(edofVec,1,24)+repmat([0 1 2 3*nely+[3 4 5 0 1 2] -3 -2 -1 3*(nelx+1)*(nely+1)+[0 1 2 3*nely+[3 4 5 0 1 2] -3 -2 -1]], nele, 1); iK = reshape(kron(edofMat,ones(24,1))',24*24*nele,1); jK = reshape(kron(edofMat,ones(1,24))',24*24*nele,1); % periodic boundary conditions n1 = [nodenrs(end, [1 end], 1) nodenrs(1, [end 1], 1) nodenrs(end, [1 end], end) nodenrs(1, [end 1], end)]; d1 = reshape([3*n1-2; 3*n1-1; 3*n1],3*numel(n1),1); n3 = [reshape(squeeze(nodenrs(end,1,2:end-1)),1,numel(squeeze(nodenrs(end,1,2:end-1))))... reshape(squeeze(nodenrs(1, 1, 2:end-1)),1,numel(squeeze(nodenrs(1, 1, 2:end-1))))... reshape(squeeze(nodenrs(end,2:end-1,1)),1,numel(squeeze(nodenrs(end,2:end-1,1))))... reshape(squeeze(nodenrs(1, 2:end-1, 1)),1,numel(squeeze(nodenrs(1, 2:end-1, 1))))... reshape(squeeze(nodenrs(2:end-1, 1, 1)),1,numel(squeeze(nodenrs(2:end-1, 1, 1))))... reshape(squeeze(nodenrs(2:end-1,1,end)),1,numel(squeeze(nodenrs(2:end-1,1,end))))... reshape(squeeze(nodenrs(2:end-1, 2:end-1, 1)),1,numel(squeeze(nodenrs(2:end-1, 2:end-1, 1))))... reshape(squeeze(nodenrs(2:end-1, 1, 2:end-1)),1,numel(squeeze(nodenrs(2:end-1, 1, 2:end-1))))... reshape(squeeze(nodenrs(end,2:end-1,2:end-1)),1,numel(squeeze(nodenrs(end,2:end-1,2:end-1))))]; d3 = reshape([3*n3-2; 3*n3-1; 3*n3],3*numel(n3),1); n4 = [reshape(squeeze(nodenrs(1, end, 2:end-1)),1,numel(squeeze(nodenrs(1, end, 2:end-1))))... reshape(squeeze(nodenrs(end,end,2:end-1)),1,numel(squeeze(nodenrs(end,end,2:end-1))))... reshape(squeeze(nodenrs(1, 2:end-1, end)),1,numel(squeeze(nodenrs(1, 2:end-1, end))))... reshape(squeeze(nodenrs(end,2:end-1,end)),1,numel(squeeze(nodenrs(end,2:end-1,end))))... reshape(squeeze(nodenrs(2:end-1,end,end)),1,numel(squeeze(nodenrs(2:end-1,end,end))))... reshape(squeeze(nodenrs(2:end-1, end, 1)),1,numel(squeeze(nodenrs(2:end-1, end, 1))))... reshape(squeeze(nodenrs(2:end-1,2:end-1,end)),1,numel(squeeze(nodenrs(2:end-1,2:end-1,end))))... reshape(squeeze(nodenrs(2:end-1,end,2:end-1)),1,numel(squeeze(nodenrs(2:end-1,end,2:end-1))))... reshape(squeeze(nodenrs(1, 2:end-1, 2:end-1)),1,numel(squeeze(nodenrs(1, 2:end-1, 2:end-1))))]; d4 = reshape([3*n4-2; 3*n4-1; 3*n4],3*numel(n4),1); n2 = setdiff(nodenrs(:),[n1(:);n3(:);n4(:)]); d2 = reshape([3*n2-2; 3*n2- 1; 3*n2],3*numel(n2),1); [cx, cy, cz] = meshgrid(0:dx:nelx, 0:dx:nely, 0:dx:nelz); % figure; scatter3(cx(:), cy(:), cz(:)); hold on; % scatter3(cx(n1(:)), cy(n1(:)), cz(n1(:)),'filled'); hold on; % scatter3(cx(n4(:)), cy(n4(:)), cz(n4(:)),'filled'); hold on; % the imposing of six linearly independent unit test strains e = eye(6); ufixed = zeros(24,6); vert_cor = [0, nelx, nelx, 0, 0, nelx, nelx, 0; 0, 0, nely, nely, 0, 0, nely, nely; 0, 0, 0, 0, nelz, nelz, nelz, nelz]; for i = 1:6 epsilon = [ e(i,1), e(i,4)/2, e(i,6)/2; e(i,4)/2, e(i,2), e(i,5)/2; e(i,6)/2, e(i,5)/2, e(i,3)]; ufixed(:,i) = reshape(epsilon*vert_cor,24,1); end % 3D boundary constraint equations wfixed = [repmat(ufixed(7:9,:),numel(squeeze(nodenrs(end,1,2:end-1))),1); repmat(ufixed( 4:6,:)-ufixed(10:12,:),numel(squeeze(nodenrs(1, 1, 2:end-1))),1); repmat(ufixed(22:24,:),numel(squeeze(nodenrs(end,2:end-1,1))),1); repmat(ufixed(13:15,:)-ufixed(10:12,:),numel(squeeze(nodenrs(1, 2:end-1, 1))),1); repmat(ufixed(16:18,:),numel(squeeze(nodenrs(2:end-1, 1, 1))),1); repmat(ufixed( 4:6,:)-ufixed(13:15,:),numel(squeeze(nodenrs(2:end-1,1,end))),1); repmat(ufixed(13:15,:),numel(squeeze(nodenrs(2:end-1, 2:end-1,1))),1); repmat(ufixed(4:6,:),numel(squeeze(nodenrs(2:end-1, 1, 2:end- 1))),1); repmat(ufixed(10:12,:),numel(squeeze(nodenrs(end,2:end-1,2:end-1))),1)]; % FE-ANALYSIS sK = reshape(KE(:)*(Emin+x(:)'.^penal*(1-Emin)),24*24*nele,1); K = sparse(iK,jK,sK); K = (K+K')/2; Kr = [K(d2,d2), K(d2,d3)+K(d2,d4); K(d3,d2)+K(d4,d2),K(d3,d3)+K(d4,d3)+K(d3,d4)+K(d4,d4)]; U(d1,:) = ufixed; U([d2;d3],:) = Kr\(-[K(d2,d1);K(d3,d1)+K(d4,d1)]*ufixed-[K(d2,d4);K(d3,d4)+K(d4,d4)]*wfixed); U(d4,:) = U(d3,:)+wfixed; ref_nodes = [cx(:),cy(:),cz(:)]; for i = 1:6 def_nodes = ref_nodes + reshape(U(:,i), 3, [])'; figure; scatter3(def_nodes(:,1), def_nodes(:,2), def_nodes(:,3)); axis equal; hold on; scatter3(def_nodes(n1(:),1), def_nodes(n1(:),2), def_nodes(n1(:),3),'filled'); hold on; % scatter3(def_nodes(n4(:),1), def_nodes(n4(:),2), def_nodes(n4(:),3),'filled'); hold on; end % homog qe = cell(6, 6); Q = zeros(6, 6); cellVolume = nele; for i = 1:6 for j = 1:6 U1 = U(:,i); U2 = U(:,j); qe{i,j} = reshape(sum((U1(edofMat)*KE).*U2(edofMat),2),nely,nelx,nelz); Q(i,j) = 1/cellVolume*sum(sum(sum((Emin+x.^penal*(1-Emin)).*qe{i,j}))); end end end function Ke = elementMatVec3D(a, b, c, DH) GN_x=[-1/sqrt(3),1/sqrt(3)]; GN_y=GN_x; GN_z=GN_x; GaussWeigh=[1,1]; Ke = zeros(24,24); L = zeros(6,9); L(1,1) = 1; L(2,5) = 1; L(3,9) = 1; L(4,2) = 1; L(4,4) = 1; L(5,6) = 1; L(5,8) = 1; L(6,3) = 1; L(6,7) = 1; for i=1:length(GN_x) for j=1:length(GN_y) for k=1:length(GN_z) x = GN_x(i);y = GN_y(j);z = GN_z(k); dNx = 1/8*[-(1-y)*(1-z), (1-y)*(1-z), (1+y)*(1-z), -(1+y)*(1-z), -(1-y)*(1+z), (1-y)*(1+z), (1+y)*(1+z) -(1+y)*(1+z)]; dNy = 1/8*[-(1-x)*(1-z) -(1+x)*(1-z), (1+x)*(1-z), (1-x)*(1-z), -(1-x)*(1+z), -(1+x)*(1+z), (1+x)*(1+z), (1-x)*(1+z)]; dNz = 1/8*[-(1-x)*(1-y) -(1+x)*(1-y) -(1+x)*(1+y) -(1-x)*(1+y), (1-x)*(1-y), (1+x)*(1-y), (1+x)*(1+y), (1-x)*(1+y)]; J = [dNx;dNy;dNz]*[ -a a a -a -a a a -a ; -b -b b b -b -b b b; -c -c -c -c c c c c]'; G = [inv(J) zeros(3) zeros(3);zeros(3) inv(J) zeros(3);zeros(3) zeros(3) inv(J)]; dN(1,1:3:24) = dNx; dN(2,1:3:24) = dNy; dN(3,1:3:24) = dNz; dN(4,2:3:24) = dNx; dN(5,2:3:24) = dNy; dN(6,2:3:24) = dNz; dN(7,3:3:24) = dNx; dN(8,3:3:24) = dNy; dN(9,3:3:24) = dNz; Be = L*G*dN; Ke = Ke + GaussWeigh(i)*GaussWeigh(j)*GaussWeigh(k)*det(J)*(Be'*DH*Be); end end end end