FEATool Multiphysics
v1.17.1
Finite Element Analysis Toolbox
|
EX_POISSON6 3D Poisson equation example on a unit cube.
[ FEA, OUT ] = EX_POISSON6( VARARGIN ) Poisson equation on a [0..1]^3 unit cube.
Input Value/{Default} Description ----------------------------------------------------------------------------------- igrid scalar 1/{0} Cell type (0=hexahedra, 1=tetrahedra) hmax scalar {1/9} Max grid cell size sfun string {sflag1} Shape function iphys scalar 0/{1} Use physics mode to define problem (=1) or directly define fea.eqn/bdr fields (=0) iplot scalar 0/{1} Plot solution (=1) . Output Value/(Size) Description ----------------------------------------------------------------------------------- fea struct Problem definition struct out struct Output struct
cOptDef = { ... 'igrid', 0; ... 'hmax', 1/9; ... 'sfun', 'sflag1'; ... 'iphys', 1; ... 'icub', 2; ... 'iplot', 1; ... 'tol', 0.05; ... 'fid', 1 }; [got,opt] = parseopt(cOptDef,varargin{:}); fid = opt.fid; % Geometry definition. gobj = gobj_block(); fea.geom.objects = { gobj }; % Grid generation. switch opt.igrid case -1 fea.grid = blockgrid(round(1/opt.hmax)); fea.grid = hex2tet(fea.grid); case 0 fea.grid = blockgrid(round(1/opt.hmax)); case 1 fea.grid = gridgen(fea,'hmax',opt.hmax,'fid',fid,'dprim',false); end n_bdr = max(fea.grid.b(3,:)); % Number of boundaries. % Problem definition. fea.sdim = { 'x' 'y' 'z' }; % Coordinate names. if ( opt.iphys==1 ) fea = addphys(fea,@poisson); % Add Poisson equation physics mode. fea.phys.poi.sfun = { opt.sfun }; % Set shape function. fea.phys.poi.eqn.coef{3,4} = { 1 }; % Set source term coefficient. fea.phys.poi.bdr.coef{1,end} = repmat({0},1,n_bdr); % Set Dirichlet boundary coefficient to zero. fea = parsephys(fea); % Check and parse physics modes. else fea.dvar = { 'u' }; % Dependent variable name. fea.sfun = { opt.sfun }; % Shape function. % Define equation system. fea.eqn.a.form = { [2 3 4;2 3 4] }; % First row indicates test function space (2=x-derivative + 3=y-derivative), % second row indicates trial function space (2=x-derivative + 3=y-derivative). fea.eqn.a.coef = { 1 }; % Coefficient used in assembling stiffness matrix. fea.eqn.f.form = { 1 }; % Test function space to evaluate in right hand side (1=function values). fea.eqn.f.coef = { 1 }; % Coefficient used in right hand side. % Define boundary conditions. fea.bdr.d = cell(1,n_bdr); [fea.bdr.d{:}] = deal(0); % Assign zero to all boundaries (Dirichlet). fea.bdr.n = cell(1,n_bdr); % No Neumann boundaries ('fea.bdr.n' empty). end % Parse and solve problem. fea = parseprob(fea); % Check and parse problem struct. fea.sol.u = solvestat(fea,'fid',fid,'icub',opt.icub); % Call to stationary solver. % Postprocessing. if ( opt.iplot>0 ) p = l_xgrid(fea.grid,opt.sfun); u_ref = refsol_poi3dcube(p(:,1),p(:,2),p(:,3),3); figure subplot(1,2,1) postplot(fea,'surfexpr','u','selexpr','(y>0.5)','axequal','on') title('Solution u') subplot(1,2,2) postplot(fea,u_ref(1:size(fea.grid.p,2)),'surfexpr','u','selexpr','(y>0.5)','axequal','on') title('Reference Solution') end % Error checking. x = linspace( 0+sqrt(eps), 1-sqrt(eps), 11 ); [x,y,z] = ndgrid(x,x,x); u = evalexpr( 'u', [x(:) y(:) z(:)]', fea ); u_ref = refsol_poi3dcube( x(:), y(:), z(:), 6 ); err = sqrt(sum((u-u_ref).^2)/sum(u_ref.^2)); if( ~isempty(fid) ) fprintf(fid,'\nL2 Error: %f\n',err) fprintf(fid,'\n\n') end out.err = err; out.pass = out.err<opt.tol; if ( nargout==0 ) clear fea out end %-----------------------------------------------