rbmatlab 0.10.01
|
00001 function U = fem_laplace(grid,params) 00002 %function U = fem_laplace(grid,params) 00003 % 00004 % function solving the elliptic laplace problem with general dirichlet and 00005 % neuman boundary data using a nodal basis. currently cartesian 00006 % grid (rectgrid) is assumed and bilinear lagrange-elements. More general 00007 % ansatz-functions and quadratures should be implemented in the future. 00008 % 00009 % Current problem: discrete gradient of the solution is not divergence 00010 % free so must be divergence-cleaned if used as velocity field 00011 % in other simulations 00012 00013 % Bernard Haasdonk 13.4.2006: 00014 00015 %if isempty(grid) 00016 % grid=grid_geometry(params); 00017 %end; 00018 00019 if ~isa(grid,'rectgrid') 00020 error('currently only rectgrid supported for fem!!'); 00021 end; 00022 00023 % indices of points on dirichlet boundary 00024 dir_vertex = zeros(size(grid.X)); 00025 00026 if params.verbose > 5 00027 disp('entering fe_laplace'); 00028 end; 00029 00030 % search all dirichlet edges and mark both corresponding vertices 00031 if params.verbose > 5 00032 disp('searching dirichlet edges'); 00033 end; 00034 dir_edge = grid.NBI == -1; 00035 nfaces = size(grid.NBI,2); 00036 for i = 1:nfaces; 00037 di = find(dir_edge(:,i)); 00038 dir_vertex(grid.VI(di,i)) = 1; 00039 dir_vertex(grid.VI(di,mod(i,nfaces)+1)) = 1; 00040 end; 00041 dir_vind = find(dir_vertex); 00042 00043 ndof = length(grid.X); 00044 dirichlet_fct = zeros(ndof,1); 00045 00046 % determine dirichlet boundary values 00047 dirichlet_fct(dir_vind) = dirichlet_values(grid.X(dir_vind),... 00048 grid.Y(dir_vind),params); 00049 00050 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 00051 % setup stiffness matrix as a sparse matrix 00052 % S(i,j) = 1 if i==j and x_i is a dirichlet-vertex 00053 % S(i,j) = 0 if i!=j and x_i is a dirichlet-vertex 00054 % S(i,j) = int sigma * grad phi_i * grad_phi_j if x_i,j are non-dirichlet 00055 % 00056 % so first compute complete stiffness-matrix 00057 % (then kronecker-kill later on before system solving) 00058 if params.verbose > 5 00059 disp('setting up stiffness matrix'); 00060 end; 00061 %S = sparse([0]zeros(ndof)); 00062 S = spalloc(ndof,ndof,13*ndof); 00063 00064 % in the general case: compute all local matrices simultaneously 00065 % L(i,j,k) contribution of element i, local basis functions j and k 00066 % 00067 % in our case, L is not element-dependent, so only compute L(j,k) 00068 % 00069 % for cartesian grid and constant sigma the local matrix is 00070 % identical for all elements 00071 % domain: [0, dx] x [0, dy], 00072 % counting right lower corner 1 then counterclockwise 00073 % detAinv = (dx*dy)^-1 00074 % grad phi_1 = detAinv * [dy-y, -x ]' 00075 % grad phi_2 = detAinv * [ y, x ]' 00076 % grad phi_3 = detAinv * [ -y, dx-x]' 00077 % grad phi_4 = detAinv * [y-dy, x-dx]' 00078 % => L(j,k) = detAinv * sigma 00079 % [1/3dx^2+1/3dy^2) -1/3dx^2+1/6dy^2 -1/6dx^2-1/6dy^2 1/6dx^2-1/3dy^2 00080 % XXXX 1/3dx^2+1/3dy^2 1/6dx^2-1/3dy^2 -1/6dx^2-1/6dy^2 00081 % XXXX XXXX 1/3dx^2+1/3dy^2 -1/3dx^2+1/6dy^2 00082 % xxxx XXXX XXXX 1/3(dx^2+dy^2) ] 00083 % 00084 nbasisfcts = size(grid.VI,2); 00085 sigma = 1; 00086 dx = (params.xrange(2)-params.xrange(1))/params.xnumintervals; 00087 dy = (params.yrange(2)-params.yrange(1))/params.ynumintervals; 00088 00089 Lx = [1/3 -1/3 -1/6 1/6; 00090 -1/3 1/3 1/6 -1/6; 00091 -1/6 1/6 1/3 -1/3; 00092 1/6 -1/6 -1/3 1/3]; 00093 00094 Ly = [1/3 1/6 -1/6 -1/3; 00095 1/6 1/3 -1/3 -1/6; 00096 -1/6 -1/3 1/3 1/6; 00097 -1/3 -1/6 1/6 1/3]; 00098 00099 L = (dx*dy)^(-1) * sigma * (Lx * dx^2 + Ly * dy^2); 00100 00101 % accumulate local matrices to single operator matrix 00102 % avoid large element loop, instead perform loop over pairs of 00103 % local basis functions 00104 00105 for j= 1:nbasisfcts 00106 for k= 1:nbasisfcts 00107 if params.verbose > 5 00108 disp(['(j,k)=(',num2str(j),',',num2str(k),')']); 00109 end; 00110 global_dofs_j = grid.VI(:,j); 00111 global_dofs_k = grid.VI(:,k); 00112 global_dofs_ind = sub2ind(size(S),global_dofs_j,global_dofs_k); 00113 if max(global_dofs_ind)>1.6e9 00114 % the following code results in t = 707.5058 s for 300x150 grid 00115 % cut S into column-block-pieces for ommiting large-index accesses 00116 % v = lget(S,global_dofs_ind); 00117 % v = v + L(j,k); 00118 % S = lset(S,global_dofs_ind,v); 00119 % the following code results in t = 476.4325 s for 300x150 grid 00120 % for z = 1:length(global_dofs_ind) 00121 % S(global_dofs_j(z),global_dofs_k(z)) = ... 00122 % S(global_dofs_j(z),global_dofs_k(z)) + L(j,k); 00123 % end; 00124 % the following code results in t = 62s s for 300x150 grid :-) 00125 % most time spent in solution of LGS 00126 V = sparse(global_dofs_j, global_dofs_k, ... 00127 L(j,k),size(S,1),size(S,2)); 00128 S = S + V; 00129 else % direct subvector operation is possible 00130 S(global_dofs_ind) = S(global_dofs_ind) + L(j,k); 00131 % S(global_dofs_ind) = S(global_dofs_ind) + L(:,j,k); 00132 end; 00133 end; 00134 end; 00135 00136 % alternative: Loop over elements, (which is to be omitted) 00137 % nelements = length(grid.A); 00138 % for i = 1:nelements 00139 % global_dofs_i = ... 00140 % S(global_dofs_i,global_dofs_i) = ... 00141 % S(global_dofs_i,global_dofs_i) + L(i,:,:); 00142 % end; 00143 00144 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 00145 % setup right hand side 00146 % b(j) = T1 + T2 + T3 if x_j is a non-dirichlet vertex 00147 % b(j) = 0 if x_j is a dirichlet vertex 00148 % 00149 % so first compute all nontrivial quantities, then set dir_vind to 0 00150 % 00151 %%%%%%%%% source term: 00152 % T1 = int source * phi_j (2D integral) 00153 % -> zero in our case as we have no source 00154 if params.verbose > 5 00155 disp('setting up right hand side'); 00156 end; 00157 00158 T1 = zeros(ndof,1); 00159 00160 %%%%%%%%% dirichlet term: 00161 % T2 = - int sigma * grad b_dir * grad phi_j (negative of a 2D integral) 00162 % after projection of b_dir to the discrete function space 00163 % this is exactly the stiffness-matrix multiplied with the b_dir 00164 % node values 00165 T2 = - S * dirichlet_fct; 00166 00167 %%%%%%%%% neuman term: 00168 % T3 = int_gamma_neu b_neu * phi_j (1D integral) 00169 % 00170 % goal: Assume b_neu also piecewise linear -> exact integration possible 00171 % for all edges, the entries of the two corresponding nodes are modified 00172 T3 = zeros(ndof,1); 00173 00174 % search all neuman edges and generate list of start and end points 00175 neu_edge = find(grid.NBI == -2); 00176 [i,j] = ind2sub(size(grid.NBI),neu_edge); 00177 % list of global indices of start points 00178 neu_vind1 = grid.VI(neu_edge); 00179 % list of global indices of end points 00180 j2 = mod(j,nfaces)+1; 00181 neu_edge_plus1 = sub2ind(size(grid.NBI),i,j2); 00182 neu_vind2 = grid.VI(neu_edge_plus1); 00183 00184 neu_vals1= neuman_values(grid.X(neu_vind1),... 00185 grid.Y(neu_vind1),[],[],[],params); 00186 neu_vals2= neuman_values(grid.X(neu_vind2),... 00187 grid.Y(neu_vind2),[],[],[],params); 00188 % generate lists of neuman integral values 00189 % neu_int1(neu_edgenr) is the integral contribution to start node 00190 neu_lengths = grid.EL(neu_edge); 00191 00192 % determine elementary integrals 00193 % (phi_1 is basis function of start point, phi_2 of end point) 00194 % int_edge Phi_1 Phi_1 = 1/3 |S| 00195 % int_edge Phi_1 Phi_2 = 1/6 |S| 00196 % int_edge Phi_2 Phi_2 = 1/3 |S| 00197 00198 % neu_int1 = int_edge Phi_1 (neu_val1 * Phi1 + neu_val2 * Phi2) 00199 neu_int1 = neu_lengths.* (1/3 * neu_vals1 + 1/6 * neu_vals2); 00200 neu_int2 = neu_lengths.* (1/6 * neu_vals1 + 1/3 * neu_vals2); 00201 T3(neu_vind1) = T3(neu_vind1) + neu_int1; 00202 T3(neu_vind2) = T3(neu_vind2) + neu_int2; 00203 00204 % assemble 00205 b = T1 + T2 + T3; 00206 b(dir_vind) = 0; 00207 00208 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 00209 %delete dirichlet-nodes in matrix 00210 S(dir_vind,:) = 0; 00211 S(:,dir_vind) = 0; 00212 S(dir_vind,dir_vind) = speye(length(dir_vind)); 00213 00214 % solve homogeneous problem 00215 % U = S^-1 * b; 00216 if params.verbose > 5 00217 disp('solving LGS system'); 00218 end; 00219 00220 %[U, flag] = pcg(S,b); 00221 % keyboard; 00222 % nonsymmetric solvers: 00223 % [U, flag] = bicgstab(S,b,[],1000); => 65sec 00224 % [U, flag] = cgs(S,b,[],1000); % => 47 sec 00225 % symmetric solver, non pd: 00226 [U, flag] = symmlq(S,b,[],1000); % => 31 sec 00227 % [U, flag] = minres(S,b,[],1000); % => 31 sec 00228 % symmetric solver, pd: 00229 % [U, flag] = pcg(S,b,[],1000); % => 46 sec 00230 if flag>0 00231 disp(['error in system solution, solver return flag = ', ... 00232 num2str(flag)]); 00233 keyboard; 00234 end; 00235 00236 % add dirichlet boundary values 00237 U = U + dirichlet_fct; 00238 00239 00240 00241 00242 00243 % TO BE ADJUSTED TO NEW SYNTAX 00244 %| \docupdate