class ExplicitIsotropicElectroMechanics_108(Material):
    """ Electromechanical model in terms of internal energy
                        W(C,D) = W_mn(C) + 1/2/eps_2/J (FD0*FD0)
                        W_mn(C) = u1*C:I+u2*G:I - 2*(u1+2*u2)*lnJ + lamb/2*(J-1)**2
    """
    def __init__(self, ndim, **kwargs):
        mtype = type(self).__name__
        super(ExplicitIsotropicElectroMechanics_108,
              self).__init__(mtype, ndim, **kwargs)
        # REQUIRES SEPARATELY
        self.nvar = self.ndim + 1
        self.energy_type = "internal_energy"
        self.legendre_transform = LegendreTransform()
        self.nature = "nonlinear"
        self.fields = "electro_mechanics"

        if self.ndim == 2:
            self.H_VoigtSize = 5
        elif self.ndim == 3:
            self.H_VoigtSize = 9

        # LOW LEVEL DISPATCHER
        self.has_low_level_dispatcher = True
        # self.has_low_level_dispatcher = False

    def KineticMeasures(self, F, ElectricFieldx, elem=0):
        from Florence.MaterialLibrary.LLDispatch._ExplicitIsotropicElectroMechanics_108_ import KineticMeasures
        D, stress = KineticMeasures(self, np.ascontiguousarray(F),
                                    ElectricFieldx)
        return D, stress, None

    def Hessian(self,
                StrainTensors,
                ElectricDisplacementx,
                elem=0,
                gcounter=0):

        mu1 = self.mu1
        mu2 = self.mu2
        lamb = self.lamb
        eps_2 = self.eps_2

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]

        D = ElectricDisplacementx.reshape(self.ndim, 1)
        Dx = D.reshape(self.ndim)
        DD = np.dot(D.T, D)[0, 0]

        C_mech = 2.*mu2/J*(2*einsum('ij,kl',b,b) - einsum('ik,jl',b,b) - einsum('il,jk',b,b)) +\
            2.*(mu1+2.*mu2)/J*( einsum("ik,jl",I,I)+einsum("il,jk",I,I) ) + \
            lamb*(2.*J-1.)*einsum("ij,kl",I,I) - lamb*(J-1)*( einsum("ik,jl",I,I)+einsum("il,jk",I,I) )

        C_elect = 1./eps_2*(0.5*DD*(einsum('ik,jl',I,I) + einsum('il,jk',I,I) + einsum('ij,kl',I,I) ) - \
                einsum('ij,k,l',I,Dx,Dx) - einsum('i,j,kl',Dx,Dx,I))

        self.elasticity_tensor = C_mech + C_elect

        self.coupling_tensor = 1. / eps_2 * (einsum('ik,j', I, Dx) + einsum(
            'i,jk', Dx, I) - einsum('ij,k', I, Dx))

        self.dielectric_tensor = 1. / eps_2 * I

        # TRANSFORM TENSORS TO THEIR ENTHALPY COUNTERPART
        E_Voigt, P_Voigt, C_Voigt = self.legendre_transform.InternalEnergyToEnthalpy(
            self.dielectric_tensor, self.coupling_tensor,
            self.elasticity_tensor)

        # BUILD HESSIAN
        factor = -1.
        H1 = np.concatenate((C_Voigt, factor * P_Voigt), axis=1)
        H2 = np.concatenate((factor * P_Voigt.T, E_Voigt), axis=1)
        H_Voigt = np.concatenate((H1, H2), axis=0)

        return H_Voigt

    def CauchyStress(self,
                     StrainTensors,
                     ElectricDisplacementx,
                     elem=0,
                     gcounter=0):

        mu1 = self.mu1
        mu2 = self.mu2
        lamb = self.lamb
        eps_2 = self.eps_2

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]
        D = ElectricDisplacementx.reshape(self.ndim, 1)
        DD = np.dot(D.T, D)[0, 0]

        trb = trace(b)
        if self.ndim == 2:
            trb += 1.

        simga_mech = 2.0*mu1/J*b + \
            2.0*mu2/J*(trb*b - np.dot(b,b)) -\
            2.0*(mu1+2*mu2)/J*I +\
            lamb*(J-1.)*I
        sigma_electric = 1. / eps_2 * (np.dot(D, D.T) - 0.5 * DD * I)
        sigma = simga_mech + sigma_electric

        return sigma

    def ElectricDisplacementx(self,
                              StrainTensors,
                              ElectricFieldx,
                              elem=0,
                              gcounter=0):
        D = self.legendre_transform.GetElectricDisplacement(
            self, StrainTensors, ElectricFieldx, elem, gcounter)

        # SANITY CHECK FOR IMPLICIT COMPUTATUTAION OF D
        # eps_2 = self.eps_2
        # E = ElectricFieldx.reshape(self.ndim,1)
        # D_exact = eps_2*E
        # D = D_exact

        return D

    def Permittivity(self,
                     StrainTensors,
                     ElectricDisplacementx,
                     elem=0,
                     gcounter=0):
        eps_2 = self.eps_2
        I = StrainTensors['I']
        self.dielectric_tensor = 1. / eps_2 * I
        # Negative definite inverse is needed due to Legendre transform
        return -eps_2 * I

    def ElectrostaticMeasures(self, F, ElectricFieldx, elem=0):
        from Florence.MaterialLibrary.LLDispatch._IsotropicElectroMechanics_108_ import KineticMeasures
        D, _, H_Voigt = KineticMeasures(self, np.ascontiguousarray(F),
                                        ElectricFieldx)
        # H_Voigt = np.linalg.inv(np.ascontiguousarray(H_Voigt[:,-self.ndim:,-self.ndim:]))
        H_Voigt = np.ascontiguousarray(H_Voigt[:, -self.ndim:, -self.ndim:])
        return D, None, H_Voigt
class IsotropicElectroMechanics_101(Material):
    """
                Simplest electromechanical model in terms of internal energy 
                        W(C,D0) = W_neo(C) + 1/2/eps_1 (FD0*FD0)
    """
    
    def __init__(self, ndim, **kwargs):
        mtype = type(self).__name__
        super(IsotropicElectroMechanics_101, self).__init__(mtype, ndim, **kwargs)
        # REQUIRES SEPARATELY
        self.nvar = self.ndim+1
        self.energy_type = "internal_energy"
        self.legendre_transform = LegendreTransform()
        self.nature = "nonlinear"
        self.fields = "electro_mechanics"  

        if self.ndim == 2:
            self.H_VoigtSize = 5
        elif self.ndim == 3:
            self.H_VoigtSize = 9

        # LOW LEVEL DISPATCHER
        self.has_low_level_dispatcher = True
        # self.has_low_level_dispatcher = False

    def KineticMeasures(self,F,ElectricFieldx, elem=0):
        from Florence.MaterialLibrary.LLDispatch._IsotropicElectroMechanics_101_ import KineticMeasures
        return KineticMeasures(self,np.ascontiguousarray(F), ElectricFieldx)


    def Hessian(self,StrainTensors,ElectricDisplacementx,elem=0,gcounter=0):

        mu = self.mu
        lamb = self.lamb
        eps_1 = self.eps_1

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]

        # D  = ElectricDisplacementx.reshape(self.ndim,1)
        Dx = ElectricDisplacementx.reshape(self.ndim)

        self.elasticity_tensor = lamb*(2.*J-1.)*einsum("ij,kl",I,I) + \
        (mu/J - lamb*(J-1))*( einsum("ik,jl",I,I)+einsum("il,jk",I,I) )

        self.coupling_tensor = J/eps_1*(einsum('ik,j',I,Dx) + einsum('i,jk',Dx,I))

        self.dielectric_tensor = J/eps_1*I 

        # TRANSFORM TENSORS TO THEIR ENTHALPY COUNTERPART
        E_Voigt, P_Voigt, C_Voigt = self.legendre_transform.InternalEnergyToEnthalpy(self.dielectric_tensor, 
            self.coupling_tensor, self.elasticity_tensor)


        # BUILD HESSIAN
        factor = -1.
        H1 = np.concatenate((C_Voigt,factor*P_Voigt),axis=1)
        H2 = np.concatenate((factor*P_Voigt.T,E_Voigt),axis=1)
        H_Voigt = np.concatenate((H1,H2),axis=0)

        return H_Voigt



    def CauchyStress(self,StrainTensors,ElectricDisplacementx,elem=0,gcounter=0):

        mu = self.mu
        lamb = self.lamb
        eps_1 = self.eps_1

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]
        D = ElectricDisplacementx.reshape(self.ndim,1)

        return 1.0*mu/J*(b - I) + lamb*(J-1)*I + J/eps_1*np.dot(D,D.T)
        # return 1.0*mu*(b - I)  + 1./eps_1*np.dot(D,D.T)


    def ElectricDisplacementx(self,StrainTensors,ElectricFieldx,elem=0,gcounter=0):
        # D = self.legendre_transform.GetElectricDisplacement(self, StrainTensors, ElectricFieldx, elem, gcounter)

        eps_1 = self.eps_1
        J = StrainTensors['J'][gcounter]
        E = ElectricFieldx.reshape(self.ndim,1)
        D_exact = eps_1/J*E
        # print np.linalg.norm(D - D_exact)
        return D_exact

        return D
class IsotropicElectroMechanics_201(Material):
    """Electromechanical model in terms of internal energy 
            W(C,D) = W_mn(C) + 1/2/eps_1 (FD0*FD0)
            W_mn(C) = u1*C:I+u2*G:I - 2*(u1+2*u2)*lnJ + lamb/2*(J-1)**2
    """
    def __init__(self, ndim, **kwargs):
        mtype = type(self).__name__
        super(IsotropicElectroMechanics_201,
              self).__init__(mtype, ndim, **kwargs)
        # REQUIRES SEPARATELY
        self.nvar = self.ndim + 1
        self.energy_type = "internal_energy"
        self.legendre_transform = LegendreTransform()
        self.nature = "nonlinear"
        self.fields = "electro_mechanics"

        if self.ndim == 2:
            self.H_VoigtSize = 5
        elif self.ndim == 3:
            self.H_VoigtSize = 9

        # LOW LEVEL DISPATCHER
        self.has_low_level_dispatcher = False

    def Hessian(self,
                StrainTensors,
                ElectricDisplacementx,
                elem=0,
                gcounter=0):

        mu1 = self.mu1
        mu2 = self.mu2
        lamb = self.lamb
        eps_1 = self.eps_1

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]

        D = ElectricDisplacementx.reshape(self.ndim, 1)
        Dx = D.reshape(self.ndim)
        DD = np.dot(D.T, D)[0, 0]

        self.elasticity_tensor = 2.*mu2/J*(2*einsum('ij,kl',b,b) - einsum('ik,jl',b,b) - einsum('il,jk',b,b)) +\
            2.*(mu1+2.*mu2)/J*( einsum("ik,jl",I,I)+einsum("il,jk",I,I) ) + \
            lamb*(2.*J-1.)*einsum("ij,kl",I,I) - lamb*(J-1)*( einsum("ik,jl",I,I)+einsum("il,jk",I,I) )

        self.coupling_tensor = J / eps_1 * (einsum('ik,j', I, Dx) +
                                            einsum('i,jk', Dx, I))

        self.dielectric_tensor = J / eps_1 * I

        # TRANSFORM TENSORS TO THEIR ENTHALPY COUNTERPART
        E_Voigt, P_Voigt, C_Voigt = self.legendre_transform.InternalEnergyToEnthalpy(
            self.dielectric_tensor, self.coupling_tensor,
            self.elasticity_tensor)

        # BUILD HESSIAN
        factor = -1.
        H1 = np.concatenate((C_Voigt, factor * P_Voigt), axis=1)
        H2 = np.concatenate((factor * P_Voigt.T, E_Voigt), axis=1)
        H_Voigt = np.concatenate((H1, H2), axis=0)

        return H_Voigt

    def CauchyStress(self,
                     StrainTensors,
                     ElectricDisplacementx,
                     elem=0,
                     gcounter=0):

        mu1 = self.mu1
        mu2 = self.mu2
        lamb = self.lamb
        eps_1 = self.eps_1

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]
        D = ElectricDisplacementx.reshape(self.ndim, 1)

        if self.ndim == 2:
            trb = trace(b) + 1
        else:
            trb = trace(b)

        simga_mech = 2.0*mu1/J*b + \
            2.0*mu2/J*(trb*b - np.dot(b,b)) -\
            2.0*(mu1+2*mu2)/J*I +\
            lamb*(J-1)*I
        sigma_electric = J / eps_1 * np.dot(D, D.T)
        sigma = simga_mech + sigma_electric

        return sigma

    def ElectricDisplacementx(self,
                              StrainTensors,
                              ElectricFieldx,
                              elem=0,
                              gcounter=0):
        D = self.legendre_transform.GetElectricDisplacement(
            self, StrainTensors, ElectricFieldx, elem, gcounter)

        # SANITY CHECK FOR IMPLICIT COMPUTATUTAION OF D
        # I = StrainTensors['I']
        # J = StrainTensors['J'][gcounter]
        # E = ElectricFieldx.reshape(self.ndim,1)
        # eps_1 = self.eps_1
        # D_exact = eps_1/J*E
        # # print np.linalg.norm(D - D_exact)
        # return D_exact

        return D
class IsotropicElectroMechanics_102(Material):
    """Electromechanical model in terms of internal energy 
                        W(C,D0) = W_neo(C) + 1/2/eps_1/J (FD0*FD0)
    """
    def __init__(self, ndim, **kwargs):
        mtype = type(self).__name__
        super(IsotropicElectroMechanics_102,
              self).__init__(mtype, ndim, **kwargs)
        # REQUIRES SEPARATELY
        self.nvar = self.ndim + 1
        self.energy_type = "internal_energy"
        self.legendre_transform = LegendreTransform()
        self.nature = "nonlinear"
        self.fields = "electro_mechanics"

        if self.ndim == 2:
            self.H_VoigtSize = 5
        elif self.ndim == 3:
            self.H_VoigtSize = 9

        # LOW LEVEL DISPATCHER
        self.has_low_level_dispatcher = False

    def Hessian(self,
                StrainTensors,
                ElectricDisplacementx,
                elem=0,
                gcounter=0):

        mu = self.mu
        lamb = self.lamb
        eps_1 = self.eps_1

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]

        D = ElectricDisplacementx.reshape(self.ndim, 1)
        Dx = D.reshape(self.ndim)
        DD = np.dot(D.T, D)[0, 0]

        C_mech = lamb*(2.*J-1.)*einsum("ij,kl",I,I) + \
        (mu/J - lamb*(J-1))*( einsum("ik,jl",I,I)+einsum("il,jk",I,I) )

        C_elect = 1./eps_1*(0.5*DD*(einsum('ik,jl',I,I) + einsum('il,jk',I,I) + einsum('ij,kl',I,I) ) - \
                einsum('ij,k,l',I,Dx,Dx) - einsum('i,j,kl',Dx,Dx,I))
        self.elasticity_tensor = C_mech + C_elect

        self.coupling_tensor = 1. / eps_1 * (einsum('ik,j', I, Dx) + einsum(
            'i,jk', Dx, I) - einsum('ij,k', I, Dx))

        self.dielectric_tensor = 1. / eps_1 * I

        # TRANSFORM TENSORS TO THEIR ENTHALPY COUNTERPART
        E_Voigt, P_Voigt, C_Voigt = self.legendre_transform.InternalEnergyToEnthalpy(
            self.dielectric_tensor, self.coupling_tensor,
            self.elasticity_tensor)

        # BUILD HESSIAN
        factor = -1.
        H1 = np.concatenate((C_Voigt, factor * P_Voigt), axis=1)
        H2 = np.concatenate((factor * P_Voigt.T, E_Voigt), axis=1)
        H_Voigt = np.concatenate((H1, H2), axis=0)

        return H_Voigt

    def CauchyStress(self,
                     StrainTensors,
                     ElectricDisplacementx,
                     elem=0,
                     gcounter=0):

        mu = self.mu
        lamb = self.lamb
        eps_1 = self.eps_1

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]
        D = ElectricDisplacementx.reshape(self.ndim, 1)
        DD = np.dot(D.T, D)[0, 0]

        return 1.0 * mu / J * (b - I) + lamb * (J - 1) * I + 1 / eps_1 * (
            np.dot(D, D.T) - 0.5 * DD * I)

    def ElectricDisplacementx(self,
                              StrainTensors,
                              ElectricFieldx,
                              elem=0,
                              gcounter=0):
        D = self.legendre_transform.GetElectricDisplacement(
            self, StrainTensors, ElectricFieldx, elem, gcounter)

        # SANITY CHECK
        # eps_1 = self.eps_1
        # D_exact = eps_1*ElectricFieldx.reshape(self.ndim,1)
        # print np.linalg.norm(D_exact - D)
        # return D_exact

        return D
Esempio n. 5
0
class Multi_Piezoelectric_100(Material):
    """ 
                Piezoelectric model in terms of internal energy 
                W(C,D) = W_mn(C) + 1/2/eps_1 (D0*D0) + 1/2/eps_2/J (FD0*FD0) 
                    + u3*(FD0/sqrt(u3*eps_3)+FN)*(FD0/sqrt(u3*eps_3)+FN) + u3*HN*HN - 2*sqrt(u3/eps_3)*D0*N
                W_mn(C) = u1*C:I+u2*G:I - 2*(u1+2*u2+u3)*lnJ + lamb/2*(J-1)**2 

    """
    def __init__(self, ndim, **kwargs):
        mtype = type(self).__name__
        super(Multi_Piezoelectric_100, self).__init__(mtype, ndim, **kwargs)
        # REQUIRES SEPARATELY
        self.nvar = self.ndim + 1
        self.energy_type = "internal_energy"
        self.legendre_transform = LegendreTransform()
        self.nature = "nonlinear"
        self.fields = "electro_mechanics"

        self.is_transversely_isotropic = True
        if self.ndim == 3:
            self.H_VoigtSize = 9
        else:
            self.H_VoigtSize = 5

        # LOW LEVEL DISPATCHER
        self.has_low_level_dispatcher = True
        # self.has_low_level_dispatcher = False

    def KineticMeasures(self, F, ElectricFieldx, elem=0):
        self.mu1 = self.mu1s[elem]
        self.mu2 = self.mu2s[elem]
        self.mu3 = self.mu3s[elem]
        self.lamb = self.lambs[elem]
        self.eps_1 = self.eps_1s[elem]
        self.eps_2 = self.eps_2s[elem]
        self.eps_3 = self.eps_3s[elem]

        from Florence.MaterialLibrary.LLDispatch._Piezoelectric_100_ import KineticMeasures
        return KineticMeasures(self, np.ascontiguousarray(F), ElectricFieldx,
                               self.anisotropic_orientations[elem][:, None])

    def Hessian(self,
                StrainTensors,
                ElectricDisplacementx,
                elem=0,
                gcounter=0):

        mu1 = self.mu1s[elem]
        mu2 = self.mu2s[elem]
        mu3 = self.mu3s[elem]
        lamb = self.lambs[elem]
        eps_1 = self.eps_1s[elem]
        eps_2 = self.eps_2s[elem]
        eps_3 = self.eps_3s[elem]

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]
        F = StrainTensors['F'][gcounter]
        H = J * np.linalg.inv(F).T
        N = self.anisotropic_orientations[elem][:, None]
        D = ElectricDisplacementx.reshape(self.ndim, 1)

        FN = np.dot(F, N)[:, 0]
        HN = np.dot(H, N)[:, 0]
        innerHN = einsum('i,i', HN, HN)
        outerHN = einsum('i,j', HN, HN)
        Dx = D.reshape(self.ndim)
        DD = np.dot(D.T, D)[0, 0]

        # Iso + Aniso
        C_mech = 2.*mu2/J* ( 2.0*einsum('ij,kl',b,b) - einsum('ik,jl',b,b) - einsum('il,jk',b,b) ) + \
                2.*(mu1+2*mu2+mu3)/J * ( einsum('ik,jl',I,I) + einsum('il,jk',I,I) ) + \
                lamb*(2.*J-1.)*einsum('ij,kl',I,I) - lamb*(J-1.) * ( einsum('ik,jl',I,I) + einsum('il,jk',I,I) ) - \
                4.*mu3/J * ( einsum('ij,kl',I,outerHN) + einsum('ij,kl',outerHN,I) ) + \
                2.*mu3/J*innerHN*(2.0*einsum('ij,kl',I,I) - einsum('ik,jl',I,I) - einsum('il,jk',I,I) ) + \
                2.*mu3/J * ( einsum('il,j,k',I,HN,HN) + einsum('jl,i,k',I,HN,HN) + \
                einsum('ik,j,l',I,HN,HN) + einsum('jk,i,l',I,HN,HN) )

        C_elect = 1./eps_2*(0.5*DD*(einsum('ik,jl',I,I) + einsum('il,jk',I,I) + einsum('ij,kl',I,I) ) - \
                einsum('ij,k,l',I,Dx,Dx) - einsum('i,j,kl',Dx,Dx,I))

        self.elasticity_tensor = C_mech + C_elect


        self.coupling_tensor = 1./eps_2*(einsum('ik,j',I,Dx) + einsum('i,jk',Dx,I) - einsum('ij,k',I,Dx)) + \
                2.*J*sqrt(mu3/eps_3)*(einsum('ik,j',I,Dx) + einsum('i,jk',Dx,I)) + \
                2.*sqrt(mu3/eps_3)*(einsum('ik,j',I,FN) + einsum('i,jk',FN,I))

        self.dielectric_tensor = J / eps_1 * np.linalg.inv(
            b) + 1. / eps_2 * I + 2. * J * sqrt(mu3 / eps_3) * I

        # TRANSFORM TENSORS TO THEIR ENTHALPY COUNTERPART
        E_Voigt, P_Voigt, C_Voigt = self.legendre_transform.InternalEnergyToEnthalpy(
            self.dielectric_tensor, self.coupling_tensor,
            self.elasticity_tensor)

        # BUILD HESSIAN
        factor = -1.
        H1 = np.concatenate((C_Voigt, factor * P_Voigt), axis=1)
        H2 = np.concatenate((factor * P_Voigt.T, E_Voigt), axis=1)
        H_Voigt = np.concatenate((H1, H2), axis=0)

        return H_Voigt

    def CauchyStress(self,
                     StrainTensors,
                     ElectricDisplacementx,
                     elem=0,
                     gcounter=0):

        mu1 = self.mu1s[elem]
        mu2 = self.mu2s[elem]
        mu3 = self.mu3s[elem]
        lamb = self.lambs[elem]
        eps_1 = self.eps_1s[elem]
        eps_2 = self.eps_2s[elem]
        eps_3 = self.eps_3s[elem]

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]
        F = StrainTensors['F'][gcounter]
        H = J * np.linalg.inv(F).T
        N = self.anisotropic_orientations[elem][:, None]
        FN = np.dot(F, N)[:, 0]
        HN = np.dot(H, N)[:, 0]
        outerFN = einsum('i,j', FN, FN)
        innerHN = einsum('i,i', HN, HN)
        outerHN = einsum('i,j', HN, HN)

        D = ElectricDisplacementx.reshape(self.ndim, 1)
        Dx = D.reshape(self.ndim)
        DD = np.dot(D.T, D)[0, 0]
        D0D = np.dot(D, D.T)

        if self.ndim == 3:
            trb = trace(b)
        elif self.ndim == 2:
            trb = trace(b) + 1.

        sigma_mech = 2.*mu1/J*b + \
            2.*mu2/J*(trb*b - np.dot(b,b)) - \
            2.*(mu1+2*mu2+mu3)/J*I + \
            lamb*(J-1)*I +\
            2*mu3/J*outerFN +\
            2*mu3/J*innerHN*I - 2*mu3/J*outerHN

        sigma_electric = 1./eps_2*(D0D - 0.5*DD*I) +\
            2.*J*sqrt(mu3/eps_3)*D0D + 2*sqrt(mu3/eps_3)*(einsum('i,j',Dx,FN) + einsum('i,j',FN,Dx))

        sigma = sigma_mech + sigma_electric

        return sigma

    def ElectricDisplacementx(self,
                              StrainTensors,
                              ElectricFieldx,
                              elem=0,
                              gcounter=0):

        # THE ELECTRIC FIELD NEEDS TO BE MODFIED TO TAKE CARE OF CONSTANT TERMS
        mu1 = self.mu1s[elem]
        mu2 = self.mu2s[elem]
        mu3 = self.mu3s[elem]
        lamb = self.lambs[elem]
        eps_1 = self.eps_1s[elem]
        eps_2 = self.eps_2s[elem]
        eps_3 = self.eps_3s[elem]

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]
        F = StrainTensors['F'][gcounter]
        H = J * np.linalg.inv(F).T
        N = self.anisotropic_orientations[elem][:, None]
        FN = np.dot(F, N)
        HN = np.dot(H, N)
        E = ElectricFieldx.reshape(self.ndim, 1)
        modElectricFieldx = (E - 2. * sqrt(mu3 / eps_3) * FN +
                             2. / J * sqrt(mu3 / eps_3) * HN)

        # D = self.legendre_transform.GetElectricDisplacement(self, StrainTensors, modElectricFieldx, elem, gcounter)

        # SANITY CHECK FOR IMPLICIT COMPUTATUTAION OF D
        inverse = np.linalg.inv(J / eps_1 * np.linalg.inv(b) + 1. / eps_2 * I +
                                2. * J * sqrt(mu3 / eps_3) * I)
        D_exact = np.dot(inverse, (E - 2. * sqrt(mu3 / eps_3) * FN +
                                   2. / J * sqrt(mu3 / eps_3) * HN))
        # print np.linalg.norm(D - D_exact)
        return D_exact

        return D
Esempio n. 6
0
class IsotropicElectroMechanics_100(Material):
    """Simplest electromechanical model with no coupling in terms of internal energy 
                W(C,D0) = W_neo(C) + 1/2/eps_1* (D0*D0)
    """
    def __init__(self, ndim, **kwargs):
        mtype = type(self).__name__
        super(IsotropicElectroMechanics_100,
              self).__init__(mtype, ndim, **kwargs)
        # REQUIRES SEPARATELY
        self.nvar = self.ndim + 1
        self.energy_type = "internal_energy"
        self.legendre_transform = LegendreTransform()
        self.nature = "nonlinear"
        self.fields = "electro_mechanics"

        if self.ndim == 2:
            self.H_VoigtSize = 5
        elif self.ndim == 3:
            self.H_VoigtSize = 9

        # LOW LEVEL DISPATCHER
        self.has_low_level_dispatcher = False

    def Hessian(self,
                StrainTensors,
                ElectricDisplacementx,
                elem=0,
                gcounter=0):

        mu = self.mu
        lamb = self.lamb
        eps_1 = self.eps_1

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]

        D = ElectricDisplacementx.reshape(self.ndim, 1)
        Dx = D.reshape(self.ndim)
        DD = np.dot(D.T, D)[0, 0]

        self.elasticity_tensor = lamb*(2.*J-1.)*einsum("ij,kl",I,I) + \
            (mu/J - lamb*(J-1))*( einsum("ik,jl",I,I)+einsum("il,jk",I,I) )

        self.coupling_tensor = np.zeros((self.ndim, self.ndim, self.ndim))

        self.dielectric_tensor = J / eps_1 * np.linalg.inv(b)

        # TRANSFORM TENSORS TO THEIR ENTHALPY COUNTERPART
        E_Voigt, P_Voigt, C_Voigt = self.legendre_transform.InternalEnergyToEnthalpy(
            self.dielectric_tensor, self.coupling_tensor,
            self.elasticity_tensor)
        # BUILD HESSIAN
        factor = -1.
        H1 = np.concatenate((C_Voigt, factor * P_Voigt), axis=1)
        H2 = np.concatenate((factor * P_Voigt.T, E_Voigt), axis=1)
        H_Voigt = np.concatenate((H1, H2), axis=0)

        return H_Voigt

    def CauchyStress(self,
                     StrainTensors,
                     ElectricDisplacementx,
                     elem=0,
                     gcounter=0):

        mu = self.mu
        lamb = self.lamb

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]

        return 1.0 * mu / J * (b - I) + lamb * (J - 1) * I

    def ElectricDisplacementx(self,
                              StrainTensors,
                              ElectricFieldx,
                              elem=0,
                              gcounter=0):
        D = self.legendre_transform.GetElectricDisplacement(
            self, StrainTensors, ElectricFieldx, elem, gcounter)

        # SANITY CHECK FOR IMPLICIT COMPUTATUTAION OF D
        # I = StrainTensors['I']
        # J = StrainTensors['J'][gcounter]
        # b = StrainTensors['b'][gcounter]
        # E = ElectricFieldx.reshape(self.ndim,1)
        # eps_1 = self.eps_1
        # D_exact = eps_1/J*np.dot(b,E)
        # # print np.linalg.norm(D - D_exact)
        # return D_exact
        # print D

        return D
class IsotropicElectroMechanics_107(Material):
    """Electromechanical model in terms of internal energy 
        
            W(C,D0) = W_mn(C) + 1/2/eps_1 (D0*D0) + 1/2/eps_2/J (FD0*FD0) + W_reg(C,D0)
            W_mn(C) = u1*C:I+u2*G:I - 2*(u1+2*u2+6*mue)*lnJ + lamb/2*(J-1)**2 + mue*(C:I)**2
            W_reg(C,D_0) = 2/eps_e/mue (F:F)(FD0*FD0) + 1/mu/eps_e**2*(FD0*FD0)**2
    """
    
    def __init__(self, ndim, **kwargs):
        mtype = type(self).__name__
        super(IsotropicElectroMechanics_107, self).__init__(mtype, ndim, **kwargs)
        # REQUIRES SEPARATELY
        self.nvar = self.ndim+1
        self.energy_type = "internal_energy"
        self.legendre_transform = LegendreTransform()
        self.nature = "nonlinear"
        self.fields = "electro_mechanics"

        if self.ndim == 2:
            self.H_VoigtSize = 5
        elif self.ndim == 3:
            self.H_VoigtSize = 9

        # LOW LEVEL DISPATCHER
        self.has_low_level_dispatcher = True
        # self.has_low_level_dispatcher = False

    def KineticMeasures(self,F,ElectricFieldx, elem=0):
        from Florence.MaterialLibrary.LLDispatch._IsotropicElectroMechanics_107_ import KineticMeasures
        return KineticMeasures(self,np.ascontiguousarray(F), ElectricFieldx)


    def Hessian(self,StrainTensors,ElectricDisplacementx,elem=0,gcounter=0):

        mu1 = self.mu1
        mu2 = self.mu2
        mue = self.mue
        lamb = self.lamb
        eps_1 = self.eps_1
        eps_2 = self.eps_2
        eps_e = self.eps_e

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]

        D  = ElectricDisplacementx.reshape(self.ndim,1)
        Dx = D.reshape(self.ndim)
        DD = np.dot(D.T,D)[0,0]
        D0D = np.dot(D,D.T)

        if self.ndim==2:
            trb = trace(b) + 1
        else:
            trb = trace(b)

        C_mech = 2.*mu2/J*(2*einsum('ij,kl',b,b) - einsum('ik,jl',b,b) - einsum('il,jk',b,b)) +\
            2.*(mu1+2.*mu2+6*mue)/J*( einsum("ik,jl",I,I)+einsum("il,jk",I,I) ) + \
            lamb*(2.*J-1.)*einsum("ij,kl",I,I) - lamb*(J-1)*( einsum("ik,jl",I,I)+einsum("il,jk",I,I) ) +\
            8.*mue/J*einsum('ij,kl',b,b)

        C_elect = 1./eps_2*(0.5*DD*(einsum('ik,jl',I,I) + einsum('il,jk',I,I) + einsum('ij,kl',I,I) ) - \
                einsum('ij,k,l',I,Dx,Dx) - einsum('i,j,kl',Dx,Dx,I))

        # REGULARISATION TERMS          
        C_reg = 8.*J/eps_e*( einsum('ij,k,l',b,Dx,Dx) + einsum('i,j,kl',Dx,Dx,b) ) + \
            8.*J**3/mue/eps_e**2*einsum('i,j,k,l',Dx,Dx,Dx,Dx)

        self.elasticity_tensor = C_mech + C_elect + C_reg

        self.coupling_tensor = 1./eps_2*(einsum('ik,j',I,Dx) + einsum('i,jk',Dx,I) - einsum('ij,k',I,Dx)) + \
            8*J/eps_e*einsum('ij,k',b,Dx) + 4*J/eps_e*trb*( einsum('ik,j',I,Dx) + einsum('i,jk',Dx,I) ) +\
            4.*J**3/mue/eps_e**2*DD* ( einsum('ik,j',I,Dx) + einsum('i,jk',Dx,I) ) + 8.*J**3/mue/eps_e**2*einsum('i,j,k',Dx,Dx,Dx)

        self.dielectric_tensor = J/eps_1*np.linalg.inv(b)  + 1./eps_2*I + \
            4.*J/eps_e*trb*I + \
            8.*J**3/mue/eps_e**2*(D0D + 0.5*DD*I)

        # TRANSFORM TENSORS TO THEIR ENTHALPY COUNTERPART
        E_Voigt, P_Voigt, C_Voigt = self.legendre_transform.InternalEnergyToEnthalpy(self.dielectric_tensor, 
            self.coupling_tensor, self.elasticity_tensor)
        
        # BUILD HESSIAN
        factor = -1.
        H1 = np.concatenate((C_Voigt,factor*P_Voigt),axis=1)
        H2 = np.concatenate((factor*P_Voigt.T,E_Voigt),axis=1)
        H_Voigt = np.concatenate((H1,H2),axis=0)

        return H_Voigt



    def CauchyStress(self,StrainTensors,ElectricDisplacementx,elem=0,gcounter=0):

        mu1 = self.mu1
        mu2 = self.mu2
        mue = self.mue
        lamb = self.lamb
        eps_2 = self.eps_2
        eps_e = self.eps_e

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]
        D = ElectricDisplacementx.reshape(self.ndim,1)
        DD = np.dot(D.T,D)[0,0]
        D0D = np.dot(D,D.T)

        if self.ndim==2:
            trb = trace(b) + 1
        else:
            trb = trace(b)

        sigma_mech = 2.0*mu1/J*b + \
            2.0*mu2/J*(trb*b - np.dot(b,b)) -\
            2.0*(mu1+2*mu2+6*mue)/J*I +\
            lamb*(J-1)*I +\
            4.*mue/J*trb*b
        sigma_electric = 1./eps_2*(D0D - 0.5*DD*I)
        # REGULARISATION TERMS
        sigma_reg = 4.*J/eps_e*(DD*b + trb*D0D) + \
            4.0*J**3/mue/eps_e**2*DD*D0D
        sigma = sigma_mech + sigma_electric + sigma_reg

        return sigma

    def ElectricDisplacementx(self,StrainTensors,ElectricFieldx,elem=0,gcounter=0):
        D = self.legendre_transform.GetElectricDisplacement(self, StrainTensors, ElectricFieldx, elem, gcounter)
        return D



    def Permittivity(self,StrainTensors,ElectricDisplacementx,elem=0,gcounter=0):

        mu1 = self.mu1
        mu2 = self.mu2
        mue = self.mue
        lamb = self.lamb
        eps_1 = self.eps_1
        eps_2 = self.eps_2
        eps_e = self.eps_e

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]

        D  = ElectricDisplacementx.reshape(self.ndim,1)
        Dx = D.reshape(self.ndim)
        DD = np.dot(D.T,D)[0,0]
        D0D = np.dot(D,D.T)

        if self.ndim==2:
            trb = trace(b) + 1
        else:
            trb = trace(b)

        self.dielectric_tensor = J/eps_1*np.linalg.inv(b)  + 1./eps_2*I + \
            4.*J/eps_e*trb*I + \
            8.*J**3/mue/eps_e**2*(D0D + 0.5*DD*I)

        # TRANSFORM TENSORS TO THEIR ENTHALPY COUNTERPART
        H_Voigt = -np.linalg.inv(self.dielectric_tensor)
        return H_Voigt


    def ElectrostaticMeasures(self,F,ElectricFieldx, elem=0):
        from Florence.MaterialLibrary.LLDispatch._IsotropicElectroMechanics_107_ import KineticMeasures
        D, _, H_Voigt = KineticMeasures(self,np.ascontiguousarray(F), ElectricFieldx)
        H_Voigt = np.ascontiguousarray(H_Voigt[:,-self.ndim:,-self.ndim:])
        return D, None, H_Voigt
Esempio n. 8
0
class IsotropicElectroMechanics_109(Material):
    """ Electromechanical model in terms of internal energy
                        W(C,D) = W_mn(C) + 1/2/eps_2/J (FD0*FD0)
                        W_mn(C) = u1*C:I+u2*G:I - 2*(u1+2*u2)*lnJ + lamb/2*(J-1)**2
    """

    def __init__(self, ndim, **kwargs):
        mtype = type(self).__name__
        super(IsotropicElectroMechanics_109, self).__init__(mtype, ndim, **kwargs)
        # REQUIRES SEPARATELY
        self.nvar = self.ndim+1
        self.energy_type = "internal_energy"
        self.legendre_transform = LegendreTransform()
        self.nature = "nonlinear"
        self.fields = "electro_mechanics"

        if self.ndim == 2:
            self.H_VoigtSize = 5
        elif self.ndim == 3:
            self.H_VoigtSize = 9

        # LOW LEVEL DISPATCHER
        self.has_low_level_dispatcher = True
        # self.has_low_level_dispatcher = False

    def KineticMeasures(self,F,ElectricFieldx, elem=0):
        from Florence.MaterialLibrary.LLDispatch._IsotropicElectroMechanics_109_ import KineticMeasures
        return KineticMeasures(self,np.ascontiguousarray(F), ElectricFieldx)


    def Hessian(self,StrainTensors,ElectricDisplacementx,elem=0,gcounter=0):

        mu_v = self.mu_v
        mu1 = self.mu1
        mu2 = self.mu2
        lamb = self.lamb
        eps_2 = self.eps_2

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]

        D  = ElectricDisplacementx.reshape(self.ndim,1)
        Dx = D.reshape(self.ndim)
        DD = np.dot(D.T,D)[0,0]

        C_mech = 2.*mu2/J*(2*einsum('ij,kl',b,b) - einsum('ik,jl',b,b) - einsum('il,jk',b,b)) +\
            2.*(mu1+2.*mu2)/J*( einsum("ik,jl",I,I)+einsum("il,jk",I,I) ) + \
            lamb*(2.*J-1.)*einsum("ij,kl",I,I) - lamb*(J-1)*( einsum("ik,jl",I,I)+einsum("il,jk",I,I) )

        C_elect = 1./eps_2*(0.5*DD*(einsum('ik,jl',I,I) + einsum('il,jk',I,I) + einsum('ij,kl',I,I) ) - \
                einsum('ij,k,l',I,Dx,Dx) - einsum('i,j,kl',Dx,Dx,I))

        # C_linear = self.H_Voigt
        C_linear = H_Voigt = lamb*einsum('ij,kl',I,I)+mu_v*(einsum('ik,jl',I,I)+einsum('il,jk',I,I))

        self.elasticity_tensor = C_mech + C_elect + C_linear

        self.coupling_tensor = 1./eps_2*(einsum('ik,j',I,Dx) + einsum('i,jk',Dx,I) - einsum('ij,k',I,Dx))

        self.dielectric_tensor = 1./eps_2*I

        # TRANSFORM TENSORS TO THEIR ENTHALPY COUNTERPART
        E_Voigt, P_Voigt, C_Voigt = self.legendre_transform.InternalEnergyToEnthalpy(self.dielectric_tensor,
            self.coupling_tensor, self.elasticity_tensor)

        # BUILD HESSIAN
        factor = -1.
        H1 = np.concatenate((C_Voigt,factor*P_Voigt),axis=1)
        H2 = np.concatenate((factor*P_Voigt.T,E_Voigt),axis=1)
        H_Voigt = np.concatenate((H1,H2),axis=0)

        return H_Voigt



    def CauchyStress(self,StrainTensors,ElectricDisplacementx,elem=0,gcounter=0):

        mu_v = self.mu_v
        mu1 = self.mu1
        mu2 = self.mu2
        lamb = self.lamb
        eps_2 = self.eps_2

        I = StrainTensors['I']
        J = StrainTensors['J'][gcounter]
        b = StrainTensors['b'][gcounter]
        strain = StrainTensors['strain'][gcounter]
        D = ElectricDisplacementx.reshape(self.ndim,1)
        DD = np.dot(D.T,D)[0,0]

        trb = trace(b)
        if self.ndim==2:
            trb += 1.

        simga_mech = 2.0*mu1/J*b + \
            2.0*mu2/J*(trb*b - np.dot(b,b)) -\
            2.0*(mu1+2*mu2)/J*I +\
            lamb*(J-1.)*I
        sigma_electric = 1./eps_2*(np.dot(D,D.T) - 0.5*DD*I)


        # CHECK IF THIS IS NECESSARY
        if self.ndim == 3:
            tre = trace(strain)
        elif self.ndim == 2:
            tre = trace(strain) + 1

        # USE FASTER TRACE FUNCTION
        sigma_linear = 2.*mu_v*strain + lamb*tre*I

        sigma = simga_mech + sigma_electric + sigma_linear

        return sigma

    def ElectricDisplacementx(self,StrainTensors,ElectricFieldx,elem=0,gcounter=0):
        D = self.legendre_transform.GetElectricDisplacement(self, StrainTensors, ElectricFieldx, elem, gcounter)

        # SANITY CHECK FOR IMPLICIT COMPUTATUTAION OF D
        # I = StrainTensors['I']
        # J = StrainTensors['J'][gcounter]
        # b = StrainTensors['b'][gcounter]
        # E = ElectricFieldx.reshape(self.ndim,1)
        # eps_1 = self.eps_1
        # eps_2 = self.eps_2
        # inverse = np.linalg.inv(J/eps_1*np.linalg.inv(b) + 1./eps_2*I)
        # D_exact = np.dot(inverse,E)
        # print np.linalg.norm(D - D_exact)
        # return D_exact

        return D