%************************************************************************** %******* Simulation for Hookean dumbbell ************ %* shear flow Fokker Planck * %* * %**************************** Tony Lelievre 29/05/2007 ******************** %************************************************************************** % run Couette_FP % This is for Hookean dumbbell % % The PDE on the pdf p is (without the advection term): % dp/dt = -\sum_j div_{X_j} ( (nabla_u X_j % - \frac{1}{4We} \sum_k A_{j,k} F(X_k)) p ) % + (1/(4We)) \sum_j\sum_k A_{j,k} \partial^2 p / \partial X_j \partial X_k % Change of variable: q=p/exp(-Pi) % Let us introduce the weight function \omega= exp(-Pi)/Z % The variational formulation for q is (r is the test function): % \int_X dq/dt r \omega = % - (1/(4We)) D^t : \int_X (\nabla_X q \otimes \nabla_X r) \omega + % NABLA_u^t : \int_X (X \otimes \nabla_X r) q \omega % where NABLA_u is a block diagonal matrix % with nabla_u repeated on the diagonal % and D=[ ... -Id 2Id -Id ...]. % The nondimensional equations are: % % Re . d_t u = (1- Epsilon) d_x,x u + d_x tau % u(0,x)=0 % u(t,0)=v % u(t,1)=0 % tau=(Epsilon/We) \int (X_1 X_2) q \omega % % The velocity on the boundary is actually progressively set to v % Discretisation for q: FULL TENSOR PRODUCT clear all; % This file contains some integrals of Hermite polynomials run Ortho_HD_normalise_20 % Physical parameters % Only d=2 and n=1 are implemented d=2; % dimension of the ambiant space n=1; % number of springs T=1; % Maximal time Re=0.1; Eps=0.9; We=0.5; v=1.; % Discretization % Space I_esp=100; % number of spacesteps dx=1/I_esp; mesh=0:dx:1; % Time N=100; % number of timesteps dt=T/N; % timstep l_max=2; % Maximal degree of Hermite polynomials % Discretisation for q: FULL TENSOR PRODUCT dim=(l_max+1)*(l_max+1); disp('Dimension of the Galerkin basis for Fokker Planck:');disp(dim); % To get the tensorial index as a function of the absolute index % 0 \leq nu(1) \leq l_max get_nu=@(i) [ floor((i-1)/(l_max+1)), i-1-floor((i-1)/(l_max+1))*(l_max+1)]; % To get the absolute index as a function of the tensorial index % 1 \leq i \leq dim get_i=@(nu) 1+nu(1)*(l_max+1)+nu(2); % Matrix S D1=diag(ones(1,n*d),-d);D1=D1((d+1):(n+1)*d,:); D2=diag(ones(1,n*d),d);D2=D2(1:n*d,:); S=-D1+D2; % Matrix D D=S*S'; % Here, D=2 Id % Operators disp('Computing matrices...'); % Operators for Fokker-Planck M=zeros(dim,dim); G_de_base=zeros(dim,dim); A=zeros(dim,dim); % \int_X (1/dt q^{n+1}) r \omega % since int_P_P= Id, this is only Id M=eye(dim,dim); % G = Nabla_u : \int_X ( \nabla_X r \otimes X ) q \omega % G depends on the timestep % G=nabla_u*G_de_base where nabla_u is the off-diagonal component % of the matrix \nabla u for i=1:dim % r_i for j=1:dim % q_j % +1 : to get the indices of Ortho_HD_normalise.m nu_i=get_nu(i)+1; nu_j=get_nu(j)+1; G_de_base(i,j)=int_DP_P(nu_i(1),nu_j(1))*int_P_X_P(nu_i(2),nu_j(2)); end end % A = D : \int_X ( \nabla_X q \otimes \nabla_X r) \omega % Here, D=2 Id % D(1,1) * \int \partial_{X_1}P_{i}(x) \partial_{X_1}P_{j}(x) \omega for i=1:dim % r_i for j=1:dim % q_j nu_i=get_nu(i)+1; nu_j=get_nu(j)+1; A(i,j)=A(i,j)+D(1,1)*int_DP_DP(nu_i(1),nu_j(1))*int_P_P(nu_i(2),nu_j(2)); end end % D(2,2) * \int \partial_{X_2}P_{i}(x) \partial_{X_2}P_{j}(x) \omega for i=1:dim % r_i for j=1:dim % q_j nu_i=get_nu(i)+1; nu_j=get_nu(j)+1; A(i,j)=A(i,j)+D(2,2)*int_P_P(nu_i(1),nu_j(1))*int_DP_DP(nu_i(2),nu_j(2)); end end % Computation of \int X_1 X_2 P_{i}(x) \omega % This vector is useful to compute tau int_X_X_q_1_2=zeros(dim,1); for i=1:dim nu_i=get_nu(i)+1; int_X_X_q_1_2(i)=int_X_P(nu_i(1))*int_X_P(nu_i(2)); end; % Operators for the velocity D1=diag(ones(1,I_esp-1),-1);D1=D1(2:I_esp,:);D1=[D1,zeros(I_esp-1,1)]; D2=diag(ones(1,I_esp-1));D2=[zeros(I_esp-1,1),D2,zeros(I_esp-1,1)]; D3=diag(ones(1,I_esp-1),+1);D3=D3(1:(I_esp-1),:);D3=[zeros(I_esp-1,1),D3]; % Mass matrix M_esp=(1/6)*D1+(2/3)*D2+(1/6)*D3; M_esp=M_esp.*dx; M_esp=sparse(M_esp); MM_esp=M_esp(:,2:I_esp); % Stiffness matrix A_esp=(-1)*D1+2*D2+(-1)*D3; A_esp=A_esp./dx; A_esp=sparse(A_esp); AA_esp=A_esp(:,2:I_esp); BB_esp=Re*MM_esp./dt+(1-Eps)*AA_esp; CLL=zeros(I_esp+1,1); % Vectors % initial conditions u=zeros(I_esp+1,1); % velocity is zero q=zeros(dim,I_esp); q(1,:)=ones(1,I_esp); % equilibrium at each point tau=zeros(I_esp,1); gradtau=zeros(I_esp-1,1); nabla_u=0; % Time iterations disp('Time interations'); for t=dt:dt:T, q_old=q; u_old=u; % Computation of u gradtau=tau(2:I_esp)-tau(1:(I_esp-1)); if ((t/T)<0.1) CLL(1)=v*10*(t/T); else CLL(1)=v ; end; CL=(Re*M_esp./dt+(1-Eps)*A_esp)*CLL; F=(Re*M_esp./dt)*u-CL+gradtau; u(2:I_esp)=BB_esp\F; if ((t/T)<0.1) u(1)=v*10*(t/T); else u(1)=v; end; % computation of tau for l=1:I_esp % iteration on the cells nabla_u=(u(l+1)-u(l))/dx; nabla_u_old=(u_old(l+1)-u_old(l))/dx; % computation of q(:,l) G=nabla_u*G_de_base; G_old=nabla_u_old*G_de_base; % Crank Nicholson M_n_p_1=(1/dt)*M - 0.5*(G-A/(4*We)); M_n=(1/dt)*M + 0.5*(G_old-A/(4*We)); q(:,l)=M_n_p_1\(M_n*q_old(:,l)); % implicit Euler % M_n_p_1=(1/dt)*M - (G-A/(4*We)); % M_n=(1/dt)*M; % q(:,l)=M_n_p_1\(M_n*q_old(:,l)); % explicit Euler % M_n=(1/dt)*M+(G_old-A/(4*We)); % q(:,l)=dt*(M_n*q_old(:,l)); % Computation of tau(l) % tau = \int_X ( X \otimes X ) q \omega tau(l)=(Eps/We)*(q(:,l)'*int_X_X_q_1_2); end; % Drawings plot(mesh',u,mesh',[tau;tau(I_esp)]); axis([0 1 -1 1.2]); drawnow; end; legend('velocity','stress');