Champs_bis.py 5.64 KB
Newer Older
Antoine MARTIN's avatar
Antoine MARTIN committed
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155



from dolfin import *
import numpy as np
import matplotlib
matplotlib.use('Agg')
import matplotlib.pyplot as plt
from mpl_toolkits.mplot3d import Axes3D
from matplotlib import cm
from matplotlib.ticker import LinearLocator, FormatStrFormatter


#cf 'Contenu des fichiers' à propos de 'reference'
def champs_admissibles(tabC,tab_tau,C0,C_HH,tA,res1):
    k,ll,mm,nn=tab_tau.shape
    tabChh=np.zeros((k,k,3,3))
    tabC0=np.zeros((k,k,3,3))
    for i in range(k):
        for j in range(k):
            tabC0[i,j]=C0
            tabChh[i,j]=C_HH
    mesh = UnitSquareMesh(res1, res1, "left/right")

    class LocalC(UserExpression):
        def __init__(self, Cij, mesh):
            self.Cij = Cij
            self.Nx, self.Ny, a, b = Cij.shape
            assert a==3 and b==3, "Wrong shape"
            self.mesh = mesh
            UserExpression.__init__(self)
        def eval_cell(self, value, x, ufc_cell):
            p = Cell(mesh, ufc_cell.index).midpoint()
            i, j = int(p[0]*(self.Nx)), int(p[1]*(self.Ny))
            value[:] = self.Cij[j, i,:,:].flatten()

        def value_shape(self):
            return (3, 3)

    Cvar=LocalC(tabC,mesh)
    C0var=LocalC(tabC0,mesh)
    CHH=LocalC(tabChh,mesh)
    tau_var= LocalC(tab_tau,mesh)
    AHH= LocalC(tA,mesh)

    # class used to define the periodic boundary map
    class PeriodicBoundary(SubDomain):
        def __init__(self, vertices, tolerance=DOLFIN_EPS):
            """ vertices stores the coordinates of the 4 unit cell corners"""
            SubDomain.__init__(self, tolerance)
            self.tol = tolerance
            self.vv = vertices
            self.a1 = self.vv[1,:]-self.vv[0,:] # first vector generating periodicity
            self.a2 = self.vv[3,:]-self.vv[0,:] # second vector generating periodicity
            # check if UC vertices form indeed a parallelogram
            assert np.linalg.norm(self.vv[2, :]-self.vv[3, :] - self.a1) <= self.tol
            assert np.linalg.norm(self.vv[2, :]-self.vv[1, :] - self.a2) <= self.tol

        def inside(self, x, on_boundary):
            # return True if on left or bottom boundary AND NOT on one of the
            # bottom-right or top-left vertices
            return bool((near(x[0], self.vv[0,0] + x[1]*self.a2[0]/self.vv[3,1], self.tol) or
                        near(x[1], self.vv[0,1] + x[0]*self.a1[1]/self.vv[1,0], self.tol)) and
                        (not ((near(x[0], self.vv[1,0], self.tol) and near(x[1], self.vv[1,1], self.tol)) or
                        (near(x[0], self.vv[3,0], self.tol) and near(x[1], self.vv[3,1], self.tol)))) and on_boundary)

        def map(self, x, y):
            if near(x[0], self.vv[2,0], self.tol) and near(x[1], self.vv[2,1], self.tol): # if on top-right corner
                y[0] = x[0] - (self.a1[0]+self.a2[0])
                y[1] = x[1] - (self.a1[1]+self.a2[1])
            elif near(x[0], self.vv[1,0] + x[1]*self.a2[0]/self.vv[2,1], self.tol): # if on right boundary
                y[0] = x[0] - self.a1[0]
                y[1] = x[1] - self.a1[1]
            else:   # should be on top boundary
                y[0] = x[0] - self.a2[0]
                y[1] = x[1] - self.a2[1]

    
    def reference_problem(C,tau,C0,k):

        def eps(v):
            return as_vector([v[0].dx(0),
                            v[1].dx(1),
                            (v[0].dx(1)+v[1].dx(0))/sqrt(2)])
   

        def sigma1(v1):
            Eps1=as_vector([1.,0.,0.])
            return dot(C0,eps(v1))+dot(tau,Eps1)
        
        def sigma2(v2):
            Eps2=as_vector([0.,1.,0.])
            return dot(C0,eps(v2))+dot(tau,Eps2)

        def sigma3(v3):
            Eps3=as_vector([0.,0.,1.])
            return dot(C0,eps(v3))+dot(tau,Eps3)

        def A_CA(e1,e2,e3):
            return as_matrix([[e1[0],e2[0],e3[0]],[e1[1],e2[1],e3[1]],[e1[2],e2[2],e3[2]]])

        def B_SA(A):
            return dot(dot(C0,A-AHH)+dot(C,AHH),inv(CHH))
        
        def t(A):
            return as_matrix([[A[0,0],A[1,0],A[2,0]],[A[0,1],A[1,1],A[2,1]],[A[0,2],A[1,2],A[2,2]]])
        
        def E_CA(A):
            return dot(t(A),dot(C,A))

        def E_SA(B):
            return dot(t(B),dot(inv(C),B))
          
        Ve = VectorElement("CG", mesh.ufl_cell(), 2)
        Re = VectorElement("R", mesh.ufl_cell(), 0)
        vertices = np.array([[0, 0], [1, 0], [1, 1], [0, 1]])
        W = FunctionSpace(mesh, MixedElement([Ve,Ve,Ve, Re,Re,Re]), constrained_domain=PeriodicBoundary(vertices, tolerance=1e-10))
        V = VectorFunctionSpace(mesh,'DG',degree=1,dim=3)

        v_1,v_2,v_3,lamb_1,lamb_2,lamb_3 = TestFunctions(W)
        dv1,dv2,dv3, dlamb1,dlamb2,dlamb3 = TrialFunctions(W)
        w = Function(W)

    
        F = inner(sigma1(dv1),eps(v_1))*dx+inner(sigma2(dv2),eps(v_2))*dx+inner(sigma3(dv3),eps(v_3))*dx
        a, L = lhs(F), rhs(F)
        a += dot(lamb_1,dv1)*dx + dot(dlamb1,v_1)*dx+dot(lamb_2,dv2)*dx+dot(dlamb2,v_2)*dx+dot(lamb_3,dv3)*dx+dot(dlamb3,v_3)*dx
        solve(a == L, w)#, solver_parameters={"linear_solver": "mumps"})#, solver_parameters={"linear_solver": "bicgstab", "preconditioner":"ilu"})
        (v1,v2,v3,lamb1,lamb2,lamb3) = split(w)              
        e1=eps(v1)+as_vector([1.,0.,0.])  
        e2=eps(v2)+as_vector([0.,1.,0.])
        e3=eps(v3)+as_vector([0.,0.,1.])
        A=A_CA(e1,e2,e3)
        B=B_SA(A)
        E1=E_CA(A)
        E2=E_SA(B)
        C_CA=np.zeros((3,3))
        S_SA=np.zeros((3,3))
        for i in range(3):
            for j in range(3):
                C_CA[i,j]=assemble(E1[i,j]*dx)
                S_SA[i,j]=assemble(E2[i,j]*dx)
        return C_CA,np.linalg.inv(S_SA)

    C_CA,C_SA = reference_problem(Cvar,tau_var,C0var,k)
    

    print(C_CA,C_SA)
    return C_CA,C_SA