You can not select more than 25 topics
			Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.
		
		
		
		
			
				
					467 lines
				
				12 KiB
			
		
		
			
		
	
	
					467 lines
				
				12 KiB
			| 
											2 years ago
										 | function [pts,tri] = xmesh2dlin(filename,varargin)
 | ||
|  | %-------------------------------------------------------------------------% | ||
|  | % XMESH2D stands for EXACT MESHER 2d. Given an arbitrary collection of | ||
|  | % NURBS curves up through polynomial degree  3 in 2D, XMESH2d will generate | ||
|  | % a higher order triangular mesh that exactly matches the input geometry. | ||
|  | % | ||
|  | % INPUT: | ||
|  | % filename: The filename of a .spline file containing the curves that | ||
|  | % define the problem geometry. | ||
|  | % | ||
|  | % Output: | ||
|  | % NODE: A nNodesx3 array containing the coordinates and weights of the | ||
|  | % nodes that make up the mesh. | ||
|  | % | ||
|  | % IEN: a 10xnel array containing the mesh connectivity information. | ||
|  | % | ||
|  | % filename.neu: xmesh2d also writes out the mesh information int he | ||
|  | % standard gambit neu trail fil e format to the file <filename>.neu | ||
|  | %-------------------------------------------------------------------------% | ||
|  | 
 | ||
|  | if nargin == 1 | ||
|  |     options.output = false; | ||
|  |     options.debug = false; | ||
|  |     options.smoothWeights = true; | ||
|  |     options.smoothNodes = false; | ||
|  |     options.thresh=1.1; | ||
|  |     options.hmax=.3; | ||
|  | 
 | ||
|  | elseif nargin == 2 | ||
|  |     options = varargin{1}; | ||
|  |     if ~exist('options.output') %#ok<*EXIST> | ||
|  |         options.output = false; | ||
|  |     end | ||
|  |     if ~exist('options.debug') | ||
|  |         options.debug = false; | ||
|  |     end | ||
|  |     if ~isfield(options,'smoothWeights') | ||
|  |         options.smoothWeights = true; | ||
|  |     end | ||
|  |     if ~exist('options.smoothNodes') | ||
|  |         options.smoothNodes = false; | ||
|  |     end | ||
|  |     if ~isfield(options,'thresh') | ||
|  |         options.thresh= 1.1; | ||
|  |     end | ||
|  |     if ~isfield(options,'hmax') | ||
|  |         options.hmax=.3; | ||
|  |     end | ||
|  | end | ||
|  | 
 | ||
|  | % Input Parser | ||
|  | [P,KV,~,~,~,face] = splineFileIn(filename); | ||
|  | 
 | ||
|  | % Normalize the knot vectors to span [0 1] | ||
|  | for kk = 1:numel(KV) | ||
|  |     KV{kk} = KV{kk}/KV{kk}(end); | ||
|  | end | ||
|  | 
 | ||
|  | % Subdividing the inputted NURBS curves into polygons | ||
|  | [node,edge,kvloc] = NURBS2poly(P,KV,options.thresh); | ||
|  | hdata.hmax=options.hmax*(max(node(:,1))-min(node(:,1))); | ||
|  | % Feeding the polygons into mesh2d | ||
|  | [pts,tri,~] = xmeshfaces(node,edge,P,KV,kvloc,face,hdata,options); | ||
|  | return | ||
|  | 
 | ||
|  | function [node,linnode,BFLAG,CFLAG]= elevateMesh(pts,tri,bNode,bflag,bc,p)
 | ||
|  | 
 | ||
|  | side = [1 2; 2 3; 3 1]; | ||
|  | side10 = [1 4 5 2; 2 6 7 3; 3 8 9 1]; | ||
|  | node = cell(length(tri),1); | ||
|  | linnode = cell(length(tri),1); | ||
|  | 
 | ||
|  | BFLAG = zeros(numel(bNode),4); | ||
|  | ctr = 1; | ||
|  | d = 14; | ||
|  | 
 | ||
|  | % Define CFLAG as a flag if the current element is curved. | ||
|  | CFLAG = false(length(tri),1); | ||
|  | cctr=0; | ||
|  | for ee = 1:length(tri) | ||
|  |     vert = pts(tri(ee,:),:); | ||
|  |     linnode{ee} = round(gen_net(vert)*10^d)/10^d; | ||
|  |     node{ee} = round(gen_net(vert)*10^d)/10^d; | ||
|  |     %Check to see if the current triangle is a boundary triangle. | ||
|  |     for bb = 1:numel(bNode) | ||
|  |          | ||
|  |         % Check to see if the triangle has a side on the boundary. | ||
|  |         for ss = 1:3 | ||
|  |             if all(single(node{ee}(side(ss,1),1:2)) == single(bNode{bb}(1,1:2))) && ... | ||
|  |                     all(single(node{ee}(side(ss,2),1:2)) == single(bNode{bb}(4,1:2))); | ||
|  |                 node{ee}(side10(ss,:),:) = bNode{bb}; | ||
|  |                 linnode{ee}(side10(ss,:),3) = bNode{bb}(:,3); | ||
|  |                  | ||
|  |                 BFLAG(ctr,1) = ee; | ||
|  |                 BFLAG(ctr,2) = ss; | ||
|  |                 BFLAG(ctr,3) = bflag(bb,1); | ||
|  |                 BFLAG(ctr,4) = bc(bflag(bb,1),2); | ||
|  |                  | ||
|  |                 if p(bflag(bb,2))>1 | ||
|  |                     CFLAG(ee) = true; | ||
|  |                     ss=evalNURBS(node{ee}(side10(ss,:),:),[0 0 0 0 1 1 1 1],0:.1:1); | ||
|  |                     if abs(ss(:,1).^2+ss(:,2).^2-1)>1e-10 | ||
|  |                         cctr=cctr+1; | ||
|  |                     end | ||
|  |                 end | ||
|  |                 ctr=ctr+1; | ||
|  |             elseif all(single(node{ee}(side(ss,2),1:2)) == single(bNode{bb}(1,1:2))) && ... | ||
|  |                     all(single(node{ee}(side(ss,1),1:2)) == single(bNode{bb}(4,1:2))); | ||
|  |                 node{ee}(side10(ss,:),:) = flipud(bNode{bb}); | ||
|  |                 linnode{ee}(side10(ss,:),3) = flipud(bNode{bb}(:,3)); | ||
|  |                  | ||
|  |                 BFLAG(ctr,1) = ee; | ||
|  |                 BFLAG(ctr,2) = ss; | ||
|  |                 BFLAG(ctr,3) = bflag(bb,1); | ||
|  |                 BFLAG(ctr,4) = bc(bflag(bb,1),2); | ||
|  |                  | ||
|  |                 if p(bflag(bb,2))>1 | ||
|  |                     CFLAG(ee) = true; | ||
|  |                     ss=evalNURBS(node{ee}(side10(ss,:),:),[0 0 0 0 1 1 1 1],0:.1:1); | ||
|  |                     if abs(ss(:,1).^2+ss(:,2).^2-1)>1e-10 | ||
|  |                         cctr=cctr+1; | ||
|  |                     end | ||
|  |                 end | ||
|  |                 ctr=ctr+1; | ||
|  |             end | ||
|  |         end | ||
|  |          | ||
|  |         % Also check to see if the triangle has a vertex on the boundary. | ||
|  |         % If it does, it is part of the blending layer, so it needs to be | ||
|  |         % flagged. | ||
|  | %         for vv = 1:3 | ||
|  | %             if single(node{ee}(vv,1:2)) == single(bNode{bb}(1,1:2)) || ... | ||
|  | %                     single(node{ee}(vv,1:2)) == single(bNode{bb}(4,1:2)); | ||
|  | %                 if p(bflag(bb,2))>1 | ||
|  | %                     CFLAG(ee) = true; | ||
|  | %                 end | ||
|  | %                 ctr=ctr+1; | ||
|  | %             end | ||
|  | %         end | ||
|  |          | ||
|  |     end | ||
|  | end | ||
|  | 
 | ||
|  | return | ||
|  | 
 | ||
|  | function [node] = gen_net(vert)
 | ||
|  | % ----------------------------------gen_net------------------------------------% | ||
|  | % Generates the control net for the current triangular Bezier patch. Reads in | ||
|  | % the verticies of a triangle, and linearly interpolates the triangle to make a | ||
|  | % 10 control point bezier element. | ||
|  | 
 | ||
|  | % INPUT: | ||
|  | % vert: a matrix in the form [x1 y1; x2 y2; x3 y3] containing the verticies of a | ||
|  | %  linear triangle, ordered counter clockwise from the bottom left vertex. | ||
|  | 
 | ||
|  | %     v3 | ||
|  | %     | \ | ||
|  | %     |  \ | ||
|  | %     |   \ | ||
|  | %     |    \ | ||
|  | %     v1----v2 | ||
|  | 
 | ||
|  | % OUTPUT: | ||
|  | % node: a 10x3 array containing the control point coordinates and weights. | ||
|  | % The standard node ordering is shown below. | ||
|  | % | ||
|  | %     3 | ||
|  | %     |\ | ||
|  | %     | \ | ||
|  | %     8  7 | ||
|  | %     |   \ | ||
|  | %     |    \ | ||
|  | %     9 10  6 | ||
|  | %     |      \ | ||
|  | %     |       \ | ||
|  | %     1--4--5--2 | ||
|  | % | ||
|  | %------------------------------------------------------------------------------% | ||
|  | 
 | ||
|  | 
 | ||
|  | node(:,1:2) = [ vert(1,1:2);... | ||
|  |     vert(2,1:2);... | ||
|  |     vert(3,1:2);... | ||
|  |     vert(1,1:2) + (vert(2,1:2)-vert(1,1:2))/3;... | ||
|  |     vert(1,1:2) + (vert(2,1:2)-vert(1,1:2))*2/3;... | ||
|  |     vert(2,1:2) + (vert(3,1:2)-vert(2,1:2))/3;... | ||
|  |     vert(2,1:2) + (vert(3,1:2)-vert(2,1:2))*2/3;... | ||
|  |     vert(3,1:2) + (vert(1,1:2)-vert(3,1:2))/3;... | ||
|  |     vert(3,1:2) + (vert(1,1:2)-vert(3,1:2))*2/3]; | ||
|  | 
 | ||
|  | node(10,1:2) = node(9,1:2) + (node(6,1:2)-node(9,1:2))/2; | ||
|  | if size(vert,2) ==3 | ||
|  |     node(4:10,3) = ones(7,1); | ||
|  |     node(1:3,3) = vert(1:3,3); | ||
|  | else | ||
|  |     node(:,3) = ones(10,1); | ||
|  | end | ||
|  | 
 | ||
|  | return | ||
|  | function [P,KV,p,BFLAG,BC,FACE] = splineFileIn(filename)
 | ||
|  | 
 | ||
|  | filename = [filename,'.spline']; | ||
|  | fileID = fopen(filename,'r'); | ||
|  | 
 | ||
|  | % Read in the file header | ||
|  | for i=1:2 | ||
|  |     line = fgetl(fileID); | ||
|  | end | ||
|  | 
 | ||
|  | % Find the total number of curves | ||
|  | line = fgetl(fileID); | ||
|  | dims = sscanf(line, '%f'); | ||
|  | nCurves = dims(1); | ||
|  | 
 | ||
|  | KV = cell(1,nCurves); | ||
|  | P  = cell(1,nCurves); | ||
|  | BFLAG = cell(1,nCurves); | ||
|  | p = zeros(1,nCurves); | ||
|  | 
 | ||
|  | for cc = 1:nCurves | ||
|  |     line = fgetl(fileID); %#ok<*NASGU> | ||
|  |     dims = fscanf(fileID,'%u'); | ||
|  |     nCP = dims(1); | ||
|  |     lKV = dims(2); | ||
|  |     p(cc)   = dims(3); | ||
|  |     line=fgetl(fileID); | ||
|  |      | ||
|  |     for ii = 1:nCP | ||
|  |         line = fgetl(fileID); | ||
|  |         tmpx = sscanf(line, '%lf'); | ||
|  |         P{cc}(ii,1) = tmpx(1); P{cc}(ii,2) = tmpx(2); P{cc}(ii,3) = tmpx(3); | ||
|  |     end | ||
|  |     line = fgetl(fileID); | ||
|  |     line = fgetl(fileID); | ||
|  |     KV{cc} = sscanf(line, '%lf'); | ||
|  |     KV{cc} = KV{cc}'; | ||
|  |     line = fgetl(fileID); | ||
|  |     nKspan = length(unique(KV{cc}))-1; | ||
|  |      | ||
|  |     for bb = 1:nKspan | ||
|  |         line = fgetl(fileID); | ||
|  |         BFLAG{cc}(bb,:) = sscanf(line, '%u'); | ||
|  |     end | ||
|  |      | ||
|  | end | ||
|  | line = fgetl(fileID);line = fgetl(fileID);line = fgetl(fileID); | ||
|  | NBC = sscanf(line, '%u'); | ||
|  | 
 | ||
|  | BC = zeros(NBC,2); | ||
|  | for bb = 1:NBC | ||
|  |     line = fgetl(fileID); | ||
|  |     BC(bb,:) = sscanf(line,'%u'); | ||
|  |      | ||
|  | end | ||
|  | 
 | ||
|  | line = fgetl(fileID); | ||
|  | 
 | ||
|  | if line ~= -1 | ||
|  |     fgetl(fileID); | ||
|  |      | ||
|  |     line = fgetl(fileID); | ||
|  |     nFace = sscanf(line,'%u'); | ||
|  |     FACE = cell(1,nFace); | ||
|  |     for ff = 1:nFace | ||
|  |         fgetl(fileID); | ||
|  |         line = fgetl(fileID); | ||
|  |         nCurve = sscanf(line,'%u'); | ||
|  |         for cc = 1:nCurve | ||
|  |             line = fgetl(fileID); | ||
|  |             FACE{ff}(cc) = sscanf(line,'%u'); | ||
|  |         end | ||
|  |     end | ||
|  | else | ||
|  |     FACE{1} = 1:nCurves; | ||
|  | end | ||
|  | 
 | ||
|  | fclose(fileID); | ||
|  | return | ||
|  | 
 | ||
|  | function [NODE, IEN,BFLAG,CFLAG]= elevateMeshAndSmooth(pts,tri,bNode)
 | ||
|  | 
 | ||
|  | side = [1 2; 2 3; 3 1]; | ||
|  | side10 = [1 4 5 2; 2 6 7 3; 3 8 9 1]; | ||
|  | node = cell(length(tri),1); | ||
|  | 
 | ||
|  | ctr = 1; | ||
|  | d = 14; | ||
|  | 
 | ||
|  | 
 | ||
|  | 
 | ||
|  | 
 | ||
|  | % Loop through elements in the mesh and elevate locally. | ||
|  | for ee = 1:length(tri) | ||
|  |     vert = pts(tri(ee,:),:); | ||
|  |     node{ee} = round(gen_net(vert)*10^d)/10^d; | ||
|  |      | ||
|  | end | ||
|  | 
 | ||
|  | % Create the global NODE and IEN arrays from the local arrays. | ||
|  | addpath('~/TriGA/') | ||
|  | addpath('~/TriGA/xmesh') | ||
|  | [NODE, IEN] = gen_arrays(node); | ||
|  | 
 | ||
|  | 
 | ||
|  | bNodeList = []; | ||
|  | g = zeros(2,max(IEN(:))); | ||
|  | % Find the displacement boundary conditions. | ||
|  | for ee = 1:length(IEN) | ||
|  |     node = NODE(IEN(:,ee),:); | ||
|  |     %Check to see if the current triangle is a boundary triangle. | ||
|  |     for bb = 1:numel(bNode) | ||
|  |          | ||
|  |         % Check to see if the triangle has a side on the boundary. | ||
|  |         for ss = 1:3 | ||
|  |             if all(single(node(side(ss,1),1:2)) == single(bNode{bb}(1,1:2))) && ... | ||
|  |                     all(single(node(side(ss,2),1:2)) == single(bNode{bb}(4,1:2))); | ||
|  |                  | ||
|  |                 ctr=ctr+1; | ||
|  |                 bNodeList = [bNodeList;IEN(side10(ss,:),ee)]; | ||
|  |                 gloc = bNode{bb}(:,1:2) - node(side10(ss,:),1:2); | ||
|  |                  | ||
|  |                 g(:,IEN(side10(ss,:),ee)) = gloc'; | ||
|  |             elseif all(single(node(side(ss,2),1:2)) == single(bNode{bb}(1,1:2))) && ... | ||
|  |                     all(single(node(side(ss,1),1:2)) == single(bNode{bb}(4,1:2))); | ||
|  |                 bNodeList = [bNodeList;IEN(side10(ss,:),ee)]; | ||
|  |                 gloc = flipud(bNode{bb}(:,1:2)) - node(side10(ss,:),1:2); | ||
|  |                 g(:,IEN(side10(ss,:),ee)) = gloc'; | ||
|  | 
 | ||
|  |                 ctr=ctr+1; | ||
|  |             end | ||
|  |         end | ||
|  |          | ||
|  |     end | ||
|  | end | ||
|  | 
 | ||
|  | % Linear elasticty smoothing to make the curvilinear mesh. | ||
|  | bNodes = unique(bNodeList); | ||
|  | NODE = LE2d(NODE,IEN,bNodes,g); | ||
|  | NODE(:,3) = 1; | ||
|  | 
 | ||
|  | NODE(:,1:2) = NODE(:,1:2)+g'; | ||
|  | 
 | ||
|  | % Defining dummy BFLAG and CFLAG variables so this function will play nice  | ||
|  | % in the TriGA environment. These will need to actaully be computed if this  | ||
|  | % is found to be worthwhile. | ||
|  | BFLAG =[]; | ||
|  | CFLAG =[]; | ||
|  | return | ||
|  | 
 | ||
|  | function [newNode] = LE2d(NODE, IEN, bNodes,g)
 | ||
|  | 
 | ||
|  | % LE2d solves the 2D linear elasticity equations for the express purpose of | ||
|  | % calculating mesh deformation. Because of this, it does not account for | ||
|  | % traction and free boundary conditions. It simply claculates the global | ||
|  | % stiffness matrix, sets the forcing matric to zero and applies diriclet | ||
|  | % boundary conditions based on the geometry of the new mesh. | ||
|  | 
 | ||
|  | % Quadrature rules | ||
|  | 
 | ||
|  | [qPts, ~, w, ~]  = quadData(16); | ||
|  | nq = length(w); | ||
|  | 
 | ||
|  | % Stiffness Properties | ||
|  | E = 100*10^9; | ||
|  | nu = .3; | ||
|  | 
 | ||
|  | lambda = nu*E/((1+nu)*(1-2*nu)); | ||
|  | mu = E/(2*(1+nu)); | ||
|  | 
 | ||
|  | D = [lambda + 2*mu, lambda, 0; lambda, lambda + 2*mu, 0; 0 0 mu]; | ||
|  | 
 | ||
|  | % Mesh/ Geometry information | ||
|  | nel = length(IEN); | ||
|  | nen = 10; | ||
|  | nNodes = max(IEN(:)); | ||
|  | nsd = 2; | ||
|  | ndof = nsd*(nNodes-length(bNodes)); | ||
|  | nedof = nen*nsd; | ||
|  | 
 | ||
|  | %Initialize the global K and F matrices. | ||
|  | K = zeros(ndof); | ||
|  | F = zeros(ndof,1); | ||
|  | 
 | ||
|  | % Construct the ID array. | ||
|  | P = 1; | ||
|  | ID  = zeros(nsd,nNodes); | ||
|  | for a = 1:nNodes | ||
|  |     if find(a==bNodes) | ||
|  |         ID(1,a) = 0; | ||
|  |         ID(2,a) = 0; | ||
|  |          | ||
|  |     else | ||
|  |         ID(1,a) = P; | ||
|  |         ID(2,a) = P+1; | ||
|  |         P = P+2; | ||
|  |     end | ||
|  | end | ||
|  | 
 | ||
|  | % Generate a lookup table for the basis. | ||
|  | [N,dN_du] = evaluateBasis(qPts); | ||
|  | 
 | ||
|  | % Loop over elements. | ||
|  | for ee = 1:nel | ||
|  |      | ||
|  |     node = NODE(IEN(:,ee),:); | ||
|  |      | ||
|  |     % Set the local stiffness matrix to zero. | ||
|  |     k = zeros(nedof); | ||
|  |     % Loop over quadrature points to construct k. | ||
|  |     for qq = 1:nq | ||
|  | 
 | ||
|  |         % Construct the B matrix | ||
|  |         B = zeros(3,nedof); | ||
|  |          | ||
|  |         [~, dR_dx,~, J_det] = tri10fast(node,N(:,qq),dN_du(:,:,qq)); | ||
|  |         for a = 1:nen | ||
|  |             B(:,nsd*(a-1)+(1:2)) = [dR_dx(a,1), 0; ... | ||
|  |                 0, dR_dx(a,2);... | ||
|  |                 dR_dx(a,2), dR_dx(a,1)]; | ||
|  |         end | ||
|  |         % Add the contributions of the current quad point to the local | ||
|  |         % stiffness matrix. | ||
|  |         J_det = abs(J_det); | ||
|  |         k = k + w(qq)*B'*D*B*J_det; | ||
|  |     end | ||
|  |      | ||
|  |      | ||
|  |      | ||
|  |     % Assemble the local element stiffness matrix to the global stiffness | ||
|  |     % matrix. | ||
|  |     for a = 1:nen | ||
|  |         for b = 1:nen | ||
|  |             for i = 1:nsd | ||
|  |                 for j = 1:nsd | ||
|  |                     p = (a-1)*nsd+i; | ||
|  |                     q = (b-1)*nsd+j; | ||
|  |                     Pa  = ID(i,IEN(a,ee)); | ||
|  |                     Pb  = ID(j,IEN(b,ee)); | ||
|  |                     if Pa ~= 0 && Pb ~= 0 | ||
|  |                         K(Pa,Pb) = K(Pa,Pb) + k(p,q); | ||
|  |                          | ||
|  |                     elseif Pa~=0 && Pb ==0 | ||
|  |                         F(Pa) = F(Pa) - k(p,q)*g(j,IEN(b,ee)); | ||
|  |                     end | ||
|  |                 end | ||
|  |             end | ||
|  |         end | ||
|  |     end | ||
|  |      | ||
|  | end | ||
|  | 
 | ||
|  | % Solve the system | ||
|  | d = K\F; | ||
|  | 
 | ||
|  | newNode = zeros(size(NODE)); | ||
|  | for i = 1:nNodes | ||
|  |     xidx = ID(1,i); | ||
|  |     yidx = ID(2,i); | ||
|  |     if xidx~=0 | ||
|  |         newNode(i,1:2) = NODE(i,1:2) + [d(xidx), d(yidx)]; | ||
|  |     else | ||
|  |         newNode(i,1:2) = NODE(i,1:2); | ||
|  |     end | ||
|  | end | ||
|  | 
 | ||
|  | return |