codes copied
This commit is contained in:
14
codes/MATLAB/leo/CREATE_MESH.m
Executable file
14
codes/MATLAB/leo/CREATE_MESH.m
Executable file
@@ -0,0 +1,14 @@
|
||||
% Program to create a structured mesh using the codes of Leo Sok
|
||||
clear all; close all
|
||||
|
||||
nodes = load('LEO_files/nodes.txt');
|
||||
ux = load('LEO_files/ux.txt') ;
|
||||
uy = load('LEO_files/uy.txt') ;
|
||||
uz = load('LEO_files/uz.txt') ;
|
||||
u = sqrt(ux.^2 + uy.^2 + uz.^2);
|
||||
resol = load('LEO_files/resol.txt') ;
|
||||
dx = resol(1); dy = resol(2) ; dz = resol(3);
|
||||
|
||||
nodes_masked = maskFEM(nodes,u);
|
||||
[N,tets,faces] = meshStructTess(nodes_masked,dx,dy,dz,0,0);
|
||||
writemesh('/home/yeye/Desktop/leomesh',N,tets,faces)
|
19
codes/MATLAB/leo/maskFEM.m
Executable file
19
codes/MATLAB/leo/maskFEM.m
Executable file
@@ -0,0 +1,19 @@
|
||||
function nodes2 = maskFEM(nodes,vel)
|
||||
|
||||
a = [];
|
||||
b = [];
|
||||
c = [];
|
||||
ind = 1;
|
||||
|
||||
for i=1:length(nodes)
|
||||
if vel(i)>0
|
||||
a(ind) = nodes(i,1);
|
||||
b(ind) = nodes(i,2);
|
||||
c(ind) = nodes(i,3);
|
||||
ind = ind +1;
|
||||
end
|
||||
end
|
||||
|
||||
nodes2 = [a', b', c'];
|
||||
|
||||
|
169
codes/MATLAB/leo/meshStructTess.m
Executable file
169
codes/MATLAB/leo/meshStructTess.m
Executable file
@@ -0,0 +1,169 @@
|
||||
function [nodes, tets, faces, P] = meshStructTess(nodes, dx, dy, dz, check_mesh, plot_mesh)
|
||||
%% [nodes, tets, faces] = meshStructTess(nodes, dx, dy, dz, check_mesh, plot_mesh)
|
||||
% Generate a tessalation from a list of structured nodes.
|
||||
% input: nodes: n times 3 matrix with on the rows the coordinates of
|
||||
% the n points in the mesh
|
||||
% dx, dy, dz: the mesh-size in the directions x, y and z
|
||||
% check_mesh: if true, then it solves a Poisson problem
|
||||
% plot_mesh: if true, then it plots the mesh
|
||||
% output: nodes: m times 3 matrix with on the rows the coordinates of
|
||||
% the m <= n points in the triangulationedi
|
||||
% tets: l times 4 matrix with on the rows the tetrahedra
|
||||
% faces: k times 3 matrix with on the rows the triangles of the
|
||||
% boundary of the mesh
|
||||
% P: Transformation matrix from input nodes to output nodes.
|
||||
% Useful also for transforming node-valued functions on
|
||||
% the input nodes to node-valued functions on the output
|
||||
% nodes
|
||||
%
|
||||
% The triangulation can be plotted using tetramesh(tets,nodes)
|
||||
|
||||
|
||||
% compute the minimum and number of points in each direction
|
||||
if size(nodes,1) < 4
|
||||
error('Triangulation needs at least 4 points')
|
||||
end
|
||||
mn = min(nodes);
|
||||
xmin = mn(1);
|
||||
ymin = mn(2);
|
||||
zmin = mn(3);
|
||||
|
||||
mn = max(nodes);
|
||||
xmax = mn(1);
|
||||
ymax = mn(2);
|
||||
zmax = mn(3);
|
||||
|
||||
nx = round((xmax-xmin)/dx +1);
|
||||
ny = round((ymax-ymin)/dy +1);
|
||||
nz = round((zmax-zmin)/dz +1);
|
||||
|
||||
Nnodes = size(nodes,1);
|
||||
|
||||
|
||||
% Define tensor which consist of nodes indices, used for the creation of
|
||||
% the tetrahedra
|
||||
|
||||
nodes3d = zeros(nx,ny,nz); % preallocate
|
||||
for i=1:Nnodes
|
||||
nodes3d(round((nodes(i,1)-xmin)/dx)+1,round((nodes(i,2)-ymin)/dy)+1,round((nodes(i,3)-zmin)/dz)+1)=i;
|
||||
end
|
||||
|
||||
|
||||
disp('Creating Tetrahedra')
|
||||
|
||||
% create tetrahedral mesh in cube, which we will reuse.
|
||||
ii = 1;
|
||||
X = zeros(8,3);
|
||||
for i=0:1
|
||||
for j=0:1
|
||||
for k=0:1
|
||||
X(ii,:) = [i,j,k];
|
||||
ii = ii+1;
|
||||
end
|
||||
end
|
||||
end
|
||||
cubetet = delaunay(X);
|
||||
|
||||
% Run through the mesh
|
||||
el = 1;
|
||||
Tetrahedra = zeros(6*(nnz(nodes3d)),4); % preallocate
|
||||
|
||||
for i=1:nx-1
|
||||
for j=1:ny-1
|
||||
for k=1:nz-1
|
||||
% take [i:i+1,j:j+1,k:k+1] as cube
|
||||
nod = zeros(1,8); % perallocate
|
||||
|
||||
for l = 1:8
|
||||
% nod is vector with node indices of cube
|
||||
nod(l) = nodes3d(i + X(l,1), j + X(l,2), k + X(l,3));
|
||||
end
|
||||
|
||||
if nnz(nod) == 8 % then the cube is inside the mesh
|
||||
tet = nod(cubetet);
|
||||
else % then there is at least one point of the cube outside the mesh
|
||||
Xs = X(logical(nod),:); % take only nodes inside the mesh
|
||||
nodx = nod(logical(nod));
|
||||
if nnz(nod) == 4 % 4 nodes, check if points are coplanar
|
||||
C = cross(Xs(2,:)-Xs(1,:), Xs(3,:)-Xs(1,:));
|
||||
cop = logical(dot(C,Xs(4,:)-Xs(1,:)));
|
||||
% if cop = 0, then points are coplanar end thus no
|
||||
% tetrahedra exists.
|
||||
end
|
||||
if (nnz(nod)>4) || (nnz(nod) == 4 && cop)
|
||||
% create tetrahedra
|
||||
tet1 = delaunay(Xs);
|
||||
tet = nodx(tet1);
|
||||
else % no tetrahedra exists
|
||||
tet = [];
|
||||
end
|
||||
end
|
||||
|
||||
% add new tetrahedra to list
|
||||
Tetrahedra(el:el+size(tet,1)-1,:) = tet;
|
||||
el = el+size(tet,1);
|
||||
end
|
||||
end
|
||||
end
|
||||
|
||||
tets = Tetrahedra(1:el-1,:); % Delete extra preallocated rows.
|
||||
clear Tetrahedra
|
||||
|
||||
disp([num2str(size(tets,1)), ' tetrahedra created'])
|
||||
|
||||
% Delete nodes which are not in any tetrahedra.
|
||||
disp('Update mesh')
|
||||
contr = zeros(size(nodes,1),1);
|
||||
for i=1:size(tets,1)
|
||||
for j=1:4
|
||||
contr(tets(i,j))=1;
|
||||
end
|
||||
end
|
||||
|
||||
nodes = nodes(logical(contr),:);
|
||||
|
||||
% compute P
|
||||
P = speye(Nnodes);
|
||||
P = P(logical(contr),:);
|
||||
|
||||
disp([num2str(nnz(~contr)), ' unused nodes in triangulation deleted.'])
|
||||
|
||||
disp('Update tetrahedra')
|
||||
|
||||
% make tetrahedra compatible with new node indices
|
||||
cumcon = cumsum(~contr)';
|
||||
tets = tets - cumcon(tets);
|
||||
|
||||
% create triangles
|
||||
if size(tets,1) == 0
|
||||
warning('No tetrahedra created')
|
||||
faces = zeros(0,3);
|
||||
else
|
||||
disp('Create Triangles')
|
||||
faces = freeBoundary(triangulation(tets,nodes));
|
||||
disp([num2str(size(faces,1)), ' triangles created'])
|
||||
end
|
||||
|
||||
% checking the mesh by solving a Poisson problem
|
||||
if check_mesh
|
||||
% Builds the P1 stiffness matrix from tets and nodes
|
||||
[A,volumes]=stifness_matrixP1_3D(tets,nodes);
|
||||
% Check if element volumes may be negative
|
||||
if any(volumes<=0)
|
||||
warning('Some elements have zero or negative volume')
|
||||
end
|
||||
% solve the Poisson problem with Dirichlet BC
|
||||
A(2:end,2:end)\ones(size(A(2:end,2:end),1),1);
|
||||
disp('If there are no warnings, it probably means that the mesh is fine')
|
||||
end
|
||||
|
||||
% Plots mesh
|
||||
if plot_mesh
|
||||
tetramesh(tets,nodes)
|
||||
xlabel('x')
|
||||
ylabel('y')
|
||||
zlabel('z')
|
||||
end
|
||||
|
||||
end
|
||||
|
97
codes/MATLAB/leo/writemesh.m
Executable file
97
codes/MATLAB/leo/writemesh.m
Executable file
@@ -0,0 +1,97 @@
|
||||
function writemesh(varargin)
|
||||
%% writemesh(path, mesh)
|
||||
% Save triangulation as path.xml and path.msh
|
||||
% mesh is a struct with fields Pts, Tet, Tri
|
||||
% alernatively one can use writemesh(path, Pts, Tet, Tri)
|
||||
% Pts should by a n times 3 matrix consisting points of the mesh
|
||||
% Tet is the m times 4 matrix consisting the tetrahedra
|
||||
% Tri is the l times 3 matrix consisting the triangles at the boundary
|
||||
|
||||
if nargin > 3
|
||||
mesh.Pts=varargin{2};
|
||||
mesh.Tet=varargin{3};
|
||||
mesh.Tri=varargin{4};
|
||||
writemesh(varargin{1},mesh,varargin(nargin));
|
||||
|
||||
elseif isstruct(varargin{2})
|
||||
rootMeshFile = varargin{1};
|
||||
|
||||
% NEW FILE
|
||||
obj = [rootMeshFile,'.msh'];
|
||||
meshfile = fopen(obj,'w');
|
||||
|
||||
obj2 = [rootMeshFile,'.xml'];
|
||||
xmlfile = fopen(obj2,'w');
|
||||
|
||||
% MESH
|
||||
fprintf(meshfile,['$MeshFormat','\n']);
|
||||
fprintf(meshfile,['2.2 0 8','\n']);
|
||||
fprintf(meshfile,['$EndMeshFormat','\n']);
|
||||
|
||||
fprintf(xmlfile,['<?xml version="1.0" encoding="UTF-8"?>','\n']);
|
||||
fprintf(xmlfile,'\n');
|
||||
fprintf(xmlfile,['<dolfin xmlns:dolfin="http://www.fenicsproject.org">','\n']);
|
||||
|
||||
mesh = varargin{2};
|
||||
|
||||
Nodes = mesh.('Pts');
|
||||
mesh = rmfield(mesh,'Pts');
|
||||
|
||||
Nodes = [(1:size(Nodes,1))' Nodes(:,1:3)];
|
||||
|
||||
% POINTS
|
||||
if ~strcmp(varargin{nargin},'mute')
|
||||
disp('Write Points')
|
||||
end
|
||||
fprintf(meshfile,['$Nodes','\n']);
|
||||
fprintf(meshfile,['%i','\n'],size(Nodes,1));
|
||||
fprintf(xmlfile,[' <mesh celltype="tetrahedron" dim="3">','\n']);
|
||||
fprintf(xmlfile,[' <vertices size="%i">','\n'],size(Nodes,1));
|
||||
|
||||
|
||||
fprintf(meshfile,'%i %13.6f %13.6f %13.6f\n',Nodes');
|
||||
|
||||
Nodes(:,1) = Nodes(:,1) - 1;
|
||||
|
||||
fprintf(xmlfile,' <vertex index="%i" x="%0.16e" y="%0.16e" z="%0.16e"/>\n',Nodes');
|
||||
|
||||
fprintf(meshfile,['$EndNodes','\n']);
|
||||
fprintf(meshfile,['$Elements','\n']);
|
||||
fprintf(meshfile,['%i','\n'],size(mesh.Tet,1)+size(mesh.Tri,1));
|
||||
fprintf(xmlfile,[' </vertices>','\n']);
|
||||
fprintf(xmlfile,[' <cells size="%i">','\n'],size(mesh.Tet,1));
|
||||
|
||||
% Triangles
|
||||
|
||||
if ~strcmp(varargin{nargin},'mute')
|
||||
disp('Write Triangles')
|
||||
end
|
||||
|
||||
tri = mesh.('Tri');
|
||||
tri = [(1:size(tri,1))' 2*ones(size(tri,1),1) 2*ones(size(tri,1),1) zeros(size(tri,1),1) 2*ones(size(tri,1),1) tri(:,1:3)];
|
||||
fprintf(meshfile,'%i %i %i %i %i %i %i %i\n',tri');
|
||||
|
||||
|
||||
|
||||
% Tetrahedra
|
||||
if ~strcmp(varargin{nargin},'mute')
|
||||
disp('Write Tetrahedra')
|
||||
end
|
||||
|
||||
tet = mesh.('Tet');
|
||||
tet = [(size(tri,1)+1:size(tri,1)+size(tet,1))' 4*ones(size(tet,1),1) 2*ones(size(tet,1),1) zeros(size(tet,1),1) ones(size(tet,1),1) tet(:,1:4)];
|
||||
fprintf(meshfile,'%i %i %i %i %i %i %i %i %i\n',tet');
|
||||
|
||||
tet = mesh.('Tet');
|
||||
tet = [(0:size(tet,1)-1)' (tet(:,1:4)-1)];
|
||||
fprintf(xmlfile,' <tetrahedron index="%i" v0="%i" v1="%i" v2="%i" v3="%i"/>\n',tet');
|
||||
|
||||
|
||||
|
||||
fprintf(meshfile,['$EndElements','\n']);
|
||||
fprintf(xmlfile,' </cells>\n </mesh>\n</dolfin>\n');
|
||||
|
||||
fclose('all');
|
||||
end
|
||||
|
||||
|
Reference in New Issue
Block a user