├── README.md └── standardFFT-linear-elastic.py /README.md: -------------------------------------------------------------------------------- 1 | # FFT-based Homogenization 2 | 3 | ![githubFFT](https://user-images.githubusercontent.com/52761240/61006553-e7a93e00-a36a-11e9-8db7-56725e06e296.png) 4 | 5 | This code is an implementation of the standard FFT homogenization method developed by Moulinec and Suquet (1998) for microstructures with linear elastic isotropic phases. 6 | It computes the effective (/homogenized) stiffness of a microstructure with two phases (in the code the inclusion has the shape of a cylinder). 7 | 8 | By now, many different FFT homogenization methods have been developded, see for example 9 | 10 | [1] J. Zeman, T.W.J. de Geus, J. Vondřejc, R.H.J. Peerlings, M.G.D. Geers. A finite element perspective on non-linear FFT-based micromechanical simulations. 11 | International Journal for Numerical Methods in Engineering, Accepted, 2016. [arXiv: 1601.05970](https://arxiv.org/abs/1601.05970). 12 | [2] T.W.J. de Geus, J. Vondřejc, J. Zeman, R.H.J. Peerlings, M.G.D. Geers. Finite strain FFT-based non-linear solvers made simple. Submitted, 2016. arXiv: [1603.08893](https://arxiv.org/abs/1603.08893). 13 | 14 | See also: https://github.com/tdegeus/GooseFFT where the basis of the codes in this repositoy can be found. 15 | -------------------------------------------------------------------------------- /standardFFT-linear-elastic.py: -------------------------------------------------------------------------------- 1 | import numpy as np 2 | import scipy.sparse.linalg as sp 3 | import itertools 4 | from mayavi import mlab #scientific visualization library 5 | 6 | ################# STANDARD FFT-HOMOGENIZATION by Moulinec & Suquet (1994) ############### 7 | 8 | # ----------------------------------- GRID ------------------------------------ 9 | 10 | ndim = 3 # number of dimensions 11 | N = 31 # number of voxels (assumed equal for all directions, needs to be an odd number) 12 | ndof = ndim**2*N**3 # number of degrees-of-freedom 13 | 14 | # ---------------------- PROJECTION, TENSORS, OPERATIONS ---------------------- 15 | 16 | # tensor operations/products: np.einsum enables index notation, avoiding loops 17 | trans2 = lambda A2 : np.einsum('ijxyz ->jixyz ' ,A2) 18 | ddot42 = lambda A4,B2: np.einsum('ijklxyz,lkxyz ->ijxyz ',A4,B2) 19 | dot22 = lambda A2,B2: np.einsum('ijxyz ,jkxyz ->ikxyz ',A2,B2) 20 | dyad22 = lambda A2,B2: np.einsum('ijxyz ,klxyz ->ijklxyz',A2,B2) 21 | 22 | # identity tensor [single tensor] 23 | i = np.eye(ndim) 24 | # identity tensors [grid of tensors] 25 | I = np.einsum('ij,xyz' , i ,np.ones([N,N,N])) 26 | I4 = np.einsum('ijkl,xyz->ijklxyz',np.einsum('il,jk',i,i),np.ones([N,N,N])) 27 | I4rt = np.einsum('ijkl,xyz->ijklxyz',np.einsum('ik,jl',i,i),np.ones([N,N,N])) 28 | I4s = (I4+I4rt)/2. # symm. 4th order tensor 29 | II = dyad22(I,I) # dyadic product of 2nd order unit tensors 30 | 31 | # ------------------- PROBLEM DEFINITION / CONSTITIVE MODEL ------------------- 32 | 33 | # phase indicator: cubical inclusion of volume fraction 34 | # here: inclusion has cylinder form 35 | phase = np.zeros([N,N,N]) 36 | r = np.sqrt((0.2*N**2)/np.pi) #radius of cylinder (20% volume fraction) 37 | for i in range(N): 38 | for j in range(N): 39 | for k in range(N): 40 | if ((i-int(N/2))**2 + (k-int(N/2))**2) <= r: 41 | phase[i,j,k]=1. 42 | 43 | ## Visualization with Mayavi 44 | # X, Y, Z = np.mgrid[0:1:N*1j, 0:1:N*1j, 0:1:N*1j] 45 | # mlab.points3d(X, Y, Z, phase, color=(0.86, 0.24, 0.22), scale_factor=0.1, mode='cube') 46 | # mlab.outline(color=(0.24, 0.56, 0.71), line_width=2.7) 47 | # mlab.show() 48 | 49 | # material parameters + function to convert to grid of scalars 50 | param = lambda M1,M2: M1*np.ones([N,N,N])*(1.-phase)+M2*np.ones([N,N,N])*phase 51 | lambda1 = 10.0 52 | lambda2 = 100.0 53 | lambdas = param(lambda1, lambda2) # Lamé constants (material1, material2) [grid of scalars] 54 | mu1 = 0.25 55 | mu2 = 2.5 56 | mu = param(mu1, mu2) # shear modulus [grid of scalars] 57 | ## stiffness tensor [grid of scalars] 58 | C4 = lambdas*II+2.*mu*I4s 59 | 60 | # ------------------------------------------------------------------------------ 61 | 62 | ## projection operator [grid of tensors] 63 | delta = lambda i,j: np.float(i==j) # Dirac delta function 64 | freq = np.arange(-(N-1)/2.,+(N+1)/2.) # coordinate axis -> freq. axis 65 | lambda0 = (lambda1 + lambda2)/2 # Lamé constant for isotropic reference material 66 | mu0 = (mu1 + mu2)/2 # shear modulus for isotropic reference material 67 | const = (lambda0 + mu0)/(mu0*(lambda0 + 2*mu0)) 68 | Greens = np.zeros([ndim,ndim,ndim,ndim,N,N,N]) # Green's function in Fourier space 69 | 70 | for k,h,i,j in itertools.product(range(ndim), repeat=4): 71 | for x,y,z in itertools.product(range(N), repeat=3): 72 | q = np.array([freq[x], freq[y], freq[z]]) # frequency vector 73 | if not q.dot(q) == 0: # zero freq. -> mean 74 | Greens[k,h,i,j,x,y,z] = (1/(4*mu0*q.dot(q))*\ 75 | (delta(k,i)*q[h]*q[j]+delta(h,i)*q[k]*q[j]+\ 76 | delta(k,j)*q[h]*q[i]+delta(h,j)*q[k]*q[i]))-\ 77 | (const*((q[i]*q[j]*q[k]*q[h])/(q.dot(q))**2)) 78 | 79 | # (inverse) Fourier transform (for each tensor component in each direction) 80 | fft = lambda x: np.fft.fftshift(np.fft.fftn (np.fft.ifftshift(x),[N,N,N])) 81 | ifft = lambda x: np.fft.fftshift(np.fft.ifftn(np.fft.ifftshift(x),[N,N,N])) 82 | 83 | # inner functions to produce IB matrix, IB = I - F^{-1} *Gamma*F*C --> eps_i+1 = IB*eps_i 84 | G = lambda x: np.real(ifft(ddot42(Greens,fft(x)))).reshape(-1) 85 | Stiff_Mat = lambda x: ddot42(C4,x.reshape(ndim,ndim,N,N,N)) 86 | G_Stiff_Mat = lambda x: G(Stiff_Mat(x)) 87 | Id = lambda x: ddot42(I4,x.reshape(ndim,ndim,N,N,N)).reshape(-1) 88 | IB = lambda x: np.add(Id(x),-1.*G_Stiff_Mat(x)) 89 | 90 | # # ----------------------------- NEWTON ITERATIONS ----------------------------- 91 | 92 | # initialize stress and strain for each 6 macroscopic strain E 93 | sig = [np.zeros([ndim,ndim,N,N,N]) for _ in range(6)] 94 | eps = [np.zeros([ndim,ndim,N,N,N]) for _ in range(6)] 95 | 96 | # set macroscopic strains, total:6 (for each direction) 97 | E = [np.zeros([ndim,ndim,N,N,N]) for _ in range(6)] 98 | 99 | E[0][0][0] = 1.0 # loading in 1,1 direction 100 | E[1][1][1] = 1.0 # loading in 2,2 direction 101 | E[2][2][2] = 1.0 # loading in 3,3 direction 102 | E[3][0][1] = 1.0 # loading in 1,2 direction 103 | E[3][1][0] = 1.0 # loading in 2,1 direction (due to symmetry) 104 | E[4][1][2] = 1.0 # loading in 2,3 direction 105 | E[4][2][1] = 1.0 # loading in 3,2 direction (due to symmetry) 106 | E[5][0][2] = 1.0 # loading in 1,3 direction 107 | E[5][2][0] = 1.0 # loading in 3,1 direction (due to symmetry) 108 | 109 | iiter = [0 for _ in range(6)] 110 | 111 | # --------------- for convergence criteria ---------------- 112 | 113 | freqMat = np.zeros([ndim, 1, N, N, N]) #[grid of scalars] 114 | for j in range(ndim): 115 | for x in range(N): 116 | for y in range(N): 117 | for z in range(N): 118 | if j==0: 119 | freqMat[j,0,x,y,z] = freq[x] 120 | elif j==1: 121 | freqMat[j,0,x,y,z] = freq[y] 122 | elif j==2: 123 | freqMat[j,0,x,y,z] = freq[z] 124 | 125 | freqMat_T = trans2(freqMat) 126 | c = int((N-1)/2) # center of frequency grid 127 | # --------------------------------------------------------- 128 | 129 | for i in range(6): 130 | sigma = np.zeros([ndim,ndim,N,N,N]) 131 | eps[i] += E[i] 132 | while True: 133 | eps[i] = IB(eps[i]) 134 | sigma = Stiff_Mat(eps[i]) 135 | 136 | # ---------------- (equilibrium-based) convergence criteria ------------------------- 137 | fou_sig = fft(sigma).reshape(ndim, ndim, N,N,N) 138 | nom = np.sqrt(np.mean(np.power(dot22(freqMat_T, fou_sig),2))) #nominator 139 | denom = np.sqrt(np.real(fou_sig[0,0,c,c,c]**2 + fou_sig[1,1,c,c,c]**2 +\ 140 | fou_sig[2,2,c,c,c]**2 + fou_sig[0,1,c,c,c]**2 +\ 141 | fou_sig[1,2,c,c,c]**2 + fou_sig[0,2,c,c,c]**2)) # denominator 142 | # --------------------------------------------------------------- 143 | 144 | if nom/denom <1.e-8 and iiter[i]>0: break 145 | iiter[i] += 1 146 | 147 | # homogenized stiffness 148 | homStiffness = np.zeros([6, 6]) 149 | 150 | # homogenization operation = 1/N Σ f_i 151 | for i in range(6): 152 | sig[i] = Stiff_Mat(eps[i]) 153 | homStiffness[0][i] = round((1.0/(N**3))*np.sum(sig[i][0][0]),4) 154 | homStiffness[1][i] = round((1.0/(N**3))*np.sum(sig[i][1][1]),4) 155 | homStiffness[2][i] = round((1.0/(N**3))*np.sum(sig[i][2][2]),4) 156 | homStiffness[3][i] = round((1.0/(N**3))*np.sum(sig[i][0][1]),4) 157 | homStiffness[4][i] = round((1.0/(N**3))*np.sum(sig[i][1][2]),4) 158 | homStiffness[5][i] = round((1.0/(N**3))*np.sum(sig[i][2][0]),4) 159 | 160 | print("Homogenized Stiffness: \n", homStiffness) 161 | --------------------------------------------------------------------------------