def triangle(md,domainname,*args): """ TRIANGLE - create model mesh using the triangle package This routine creates a model mesh using TriMesh and a domain outline, to within a certain resolution where md is a @model object, domainname is the name of an Argus domain outline file, and resolution is a characteristic length for the mesh (same unit as the domain outline unit). Riftname is an optional argument (Argus domain outline) describing rifts. Usage: md=triangle(md,domainname,resolution) or md=triangle(md,domainname, resolution, riftname) Examples: md=triangle(md,'DomainOutline.exp',1000); md=triangle(md,'DomainOutline.exp',1000,'Rifts.exp'); """ #Figure out a characteristic area. Resolution is a node oriented concept (ex a 1000m resolution node would #be made of 1000*1000 area squares). if len(args)==1: resolution=args[0] riftname='' if len(args)==2: riftname=args[0] resolution=args[1] #Check that mesh was not already run, and warn user: if md.mesh.numberofelements: choice = input('This model already has a mesh. Are you sure you want to go ahead? (y/n)') if not m.strcmp(choice,'y'): print('no meshing done ... exiting') return None area = resolution**2 #Mesh using TriMesh md.mesh=mesh2d() [md.mesh.elements,md.mesh.x,md.mesh.y,md.mesh.segments,md.mesh.segmentmarkers]=TriMesh(domainname,riftname,area) md.mesh.elements=md.mesh.elements.astype(int) md.mesh.segments=md.mesh.segments.astype(int) md.mesh.segmentmarkers=md.mesh.segmentmarkers.astype(int) #Fill in rest of fields: md.mesh.numberofelements = numpy.size(md.mesh.elements,axis=0) md.mesh.numberofvertices = numpy.size(md.mesh.x) md.mesh.vertexonboundary = numpy.zeros(md.mesh.numberofvertices,bool) md.mesh.vertexonboundary[md.mesh.segments[:,0:2]-1] = True #Now, build the connectivity tables for this mesh. [md.mesh.vertexconnectivity] = NodeConnectivity(md.mesh.elements, md.mesh.numberofvertices) [md.mesh.elementconnectivity] = ElementConnectivity(md.mesh.elements, md.mesh.vertexconnectivity) return md
def contourenvelope(md,*args): """ CONTOURENVELOPE - build a set of segments enveloping a contour .exp Usage: segments=contourenvelope(md,varargin) Example: segments=contourenvelope(md,'Stream.exp'); segments=contourenvelope(md); """ #some checks if len(args)>1: raise RuntimeError("contourenvelope error message: bad usage") if len(args)==1: flags=args[0] if isinstance(flags,str): file=flags if not os.path.exists(file): raise IOError("contourenvelope error message: file '%s' not found" % file) isfile=1 elif isinstance(flags,(bool,int,float)): #do nothing for now isfile=0 else: raise TypeError("contourenvelope error message: second argument should be a file or an elements flag") #Now, build the connectivity tables for this mesh. #Computing connectivity if numpy.size(md.mesh.vertexconnectivity,axis=0)!=md.mesh.numberofvertices and numpy.size(md.mesh.vertexconnectivity,axis=0)!=md.mesh.numberofvertices2d: [md.mesh.vertexconnectivity]=NodeConnectivity(md.mesh.elements,md.mesh.numberofvertices) if numpy.size(md.mesh.elementconnectivity,axis=0)!=md.mesh.numberofelements and numpy.size(md.mesh.elementconnectivity,axis=0)!=md.mesh.numberofelements2d: [md.mesh.elementconnectivity]=ElementConnectivity(md.mesh.elements,md.mesh.vertexconnectivity) #get nodes inside profile elementconnectivity=copy.deepcopy(md.mesh.elementconnectivity) if md.mesh.dimension()==2: elements=copy.deepcopy(md.mesh.elements) x=copy.deepcopy(md.mesh.x) y=copy.deepcopy(md.mesh.y) numberofvertices=copy.deepcopy(md.mesh.numberofvertices) numberofelements=copy.deepcopy(md.mesh.numberofelements) else: elements=copy.deepcopy(md.mesh.elements2d) x=copy.deepcopy(md.mesh.x2d) y=copy.deepcopy(md.mesh.y2d) numberofvertices=copy.deepcopy(md.mesh.numberofvertices2d) numberofelements=copy.deepcopy(md.mesh.numberofelements2d) if len(args)==1: if isfile: #get flag list of elements and nodes inside the contour nodein=ContourToMesh(elements,x,y,file,'node',1) elemin=(numpy.sum(nodein(elements),axis=1)==numpy.size(elements,axis=1)) #modify element connectivity elemout=numpy.nonzero(numpy.logical_not(elemin))[0] elementconnectivity[elemout,:]=0 elementconnectivity[numpy.nonzero(m.ismember(elementconnectivity,elemout+1))]=0 else: #get flag list of elements and nodes inside the contour nodein=numpy.zeros(numberofvertices) elemin=numpy.zeros(numberofelements) pos=numpy.nonzero(flags) elemin[pos]=1 nodein[elements[pos,:]-1]=1 #modify element connectivity elemout=numpy.nonzero(numpy.logical_not(elemin))[0] elementconnectivity[elemout,:]=0 elementconnectivity[numpy.nonzero(m.ismember(elementconnectivity,elemout+1))]=0 #Find element on boundary #First: find elements on the boundary of the domain flag=copy.deepcopy(elementconnectivity) if len(args)==1: flag[numpy.nonzero(flag)]=elemin[flag[numpy.nonzero(flag)]] elementonboundary=numpy.logical_and(numpy.prod(flag,axis=1)==0,numpy.sum(flag,axis=1)>0) #Find segments on boundary pos=numpy.nonzero(elementonboundary)[0] num_segments=numpy.size(pos) segments=numpy.zeros((num_segments*3,3),int) count=0 for el1 in pos: els2=elementconnectivity[el1,numpy.nonzero(elementconnectivity[el1,:])[0]]-1 if numpy.size(els2)>1: flag=numpy.intersect1d(numpy.intersect1d(elements[els2[0],:],elements[els2[1],:]),elements[el1,:]) nods1=elements[el1,:] nods1=numpy.delete(nods1,numpy.nonzero(nods1==flag)) segments[count,:]=[nods1[0],nods1[1],el1+1] ord1=numpy.nonzero(nods1[0]==elements[el1,:])[0][0] ord2=numpy.nonzero(nods1[1]==elements[el1,:])[0][0] #swap segment nodes if necessary if ( (ord1==0 and ord2==1) or (ord1==1 and ord2==2) or (ord1==2 and ord2==0) ): temp=segments[count,0] segments[count,0]=segments[count,1] segments[count,1]=temp segments[count,0:2]=numpy.flipud(segments[count,0:2]) count+=1 else: nods1=elements[el1,:] flag=numpy.setdiff1d(nods1,elements[els2,:]) for j in range(0,3): nods=numpy.delete(nods1,j) if numpy.any(m.ismember(flag,nods)): segments[count,:]=[nods[0],nods[1],el1+1] ord1=numpy.nonzero(nods[0]==elements[el1,:])[0][0] ord2=numpy.nonzero(nods[1]==elements[el1,:])[0][0] if ( (ord1==0 and ord2==1) or (ord1==1 and ord2==2) or (ord1==2 and ord2==0) ): temp=segments[count,0] segments[count,0]=segments[count,1] segments[count,1]=temp segments[count,0:2]=numpy.flipud(segments[count,0:2]) count+=1 segments=segments[0:count,:] return segments
def extract(md, area): # {{{ """ extract - extract a model according to an Argus contour or flag list This routine extracts a submodel from a bigger model with respect to a given contour md must be followed by the corresponding exp file or flags list It can either be a domain file (argus type, .exp extension), or an array of element flags. If user wants every element outside the domain to be extract2d, add '~' to the name of the domain file (ex: '~HO.exp') an empty string '' will be considered as an empty domain a string 'all' will be considered as the entire domain Usage: md2=extract(md,area) Examples: md2=extract(md,'Domain.exp') See also: EXTRUDE, COLLAPSE """ #copy model md1 = copy.deepcopy(md) #get elements that are inside area flag_elem = FlagElements(md1, area) if not np.any(flag_elem): raise RuntimeError("extracted model is empty") #kick out all elements with 3 dirichlets spc_elem = np.nonzero(np.logical_not(flag_elem))[0] spc_node = np.unique(md1.mesh.elements[spc_elem, :]) - 1 flag = np.ones(md1.mesh.numberofvertices) flag[spc_node] = 0 pos = np.nonzero( np.logical_not(np.sum(flag[md1.mesh.elements - 1], axis=1)))[0] flag_elem[pos] = 0 #extracted elements and nodes lists pos_elem = np.nonzero(flag_elem)[0] pos_node = np.unique(md1.mesh.elements[pos_elem, :]) - 1 #keep track of some fields numberofvertices1 = md1.mesh.numberofvertices numberofelements1 = md1.mesh.numberofelements numberofvertices2 = np.size(pos_node) numberofelements2 = np.size(pos_elem) flag_node = np.zeros(numberofvertices1) flag_node[pos_node] = 1 #Create Pelem and Pnode (transform old nodes in new nodes and same thing for the elements) Pelem = np.zeros(numberofelements1, int) Pelem[pos_elem] = np.arange(1, numberofelements2 + 1) Pnode = np.zeros(numberofvertices1, int) Pnode[pos_node] = np.arange(1, numberofvertices2 + 1) #renumber the elements (some node won't exist anymore) elements_1 = copy.deepcopy(md1.mesh.elements) elements_2 = elements_1[pos_elem, :] elements_2[:, 0] = Pnode[elements_2[:, 0] - 1] elements_2[:, 1] = Pnode[elements_2[:, 1] - 1] elements_2[:, 2] = Pnode[elements_2[:, 2] - 1] if md1.mesh.__class__.__name__ == 'mesh3dprisms': elements_2[:, 3] = Pnode[elements_2[:, 3] - 1] elements_2[:, 4] = Pnode[elements_2[:, 4] - 1] elements_2[:, 5] = Pnode[elements_2[:, 5] - 1] #OK, now create the new model! #take every field from model md2 = copy.deepcopy(md1) #automatically modify fields #loop over model fields model_fields = vars(md1) for fieldi in model_fields: #get field field = getattr(md1, fieldi) fieldsize = np.shape(field) if hasattr(field, '__dict__') and not m.ismember( fieldi, ['results'])[0]: #recursive call object_fields = vars(field) for fieldj in object_fields: #get field field = getattr(getattr(md1, fieldi), fieldj) fieldsize = np.shape(field) if len(fieldsize): #size = number of nodes * n if fieldsize[0] == numberofvertices1: setattr(getattr(md2, fieldi), fieldj, field[pos_node]) elif fieldsize[0] == numberofvertices1 + 1: setattr(getattr(md2, fieldi), fieldj, np.vstack((field[pos_node], field[-1, :]))) #size = number of elements * n elif fieldsize[0] == numberofelements1: setattr(getattr(md2, fieldi), fieldj, field[pos_elem]) else: if len(fieldsize): #size = number of nodes * n if fieldsize[0] == numberofvertices1: setattr(md2, fieldi, field[pos_node]) elif fieldsize[0] == numberofvertices1 + 1: setattr(md2, fieldi, np.hstack((field[pos_node], field[-1, :]))) #size = number of elements * n elif fieldsize[0] == numberofelements1: setattr(md2, fieldi, field[pos_elem]) #modify some specific fields #Mesh md2.mesh.numberofelements = numberofelements2 md2.mesh.numberofvertices = numberofvertices2 md2.mesh.elements = elements_2 #mesh.uppervertex mesh.lowervertex if md1.mesh.__class__.__name__ == 'mesh3dprisms': md2.mesh.uppervertex = md1.mesh.uppervertex[pos_node] pos = np.where(~np.isnan(md2.mesh.uppervertex))[0] md2.mesh.uppervertex[pos] = Pnode[ md2.mesh.uppervertex[pos].astype(int) - 1] md2.mesh.lowervertex = md1.mesh.lowervertex[pos_node] pos = np.where(~np.isnan(md2.mesh.lowervertex))[0] md2.mesh.lowervertex[pos] = Pnode[ md2.mesh.lowervertex[pos].astype(int) - 1] md2.mesh.upperelements = md1.mesh.upperelements[pos_elem] pos = np.where(~np.isnan(md2.mesh.upperelements))[0] md2.mesh.upperelements[pos] = Pelem[ md2.mesh.upperelements[pos].astype(int) - 1] md2.mesh.lowerelements = md1.mesh.lowerelements[pos_elem] pos = np.where(~np.isnan(md2.mesh.lowerelements))[0] md2.mesh.lowerelements[pos] = Pelem[ md2.mesh.lowerelements[pos].astype(int) - 1] #Initial 2d mesh if md1.mesh.__class__.__name__ == 'mesh3dprisms': flag_elem_2d = flag_elem[np.arange(0, md1.mesh.numberofelements2d)] pos_elem_2d = np.nonzero(flag_elem_2d)[0] flag_node_2d = flag_node[np.arange(0, md1.mesh.numberofvertices2d)] pos_node_2d = np.nonzero(flag_node_2d)[0] md2.mesh.numberofelements2d = np.size(pos_elem_2d) md2.mesh.numberofvertices2d = np.size(pos_node_2d) md2.mesh.elements2d = md1.mesh.elements2d[pos_elem_2d, :] md2.mesh.elements2d[:, 0] = Pnode[md2.mesh.elements2d[:, 0] - 1] md2.mesh.elements2d[:, 1] = Pnode[md2.mesh.elements2d[:, 1] - 1] md2.mesh.elements2d[:, 2] = Pnode[md2.mesh.elements2d[:, 2] - 1] md2.mesh.x2d = md1.mesh.x[pos_node_2d] md2.mesh.y2d = md1.mesh.y[pos_node_2d] #Edges if m.strcmp(md.mesh.domaintype(), '2Dhorizontal'): if np.ndim(md2.mesh.edges) > 1 and np.size( md2.mesh.edges, axis=1 ) > 1: #do not use ~isnan because there are some np.nans... #renumber first two columns pos = np.nonzero(md2.mesh.edges[:, 3] != -1)[0] md2.mesh.edges[:, 0] = Pnode[md2.mesh.edges[:, 0] - 1] md2.mesh.edges[:, 1] = Pnode[md2.mesh.edges[:, 1] - 1] md2.mesh.edges[:, 2] = Pelem[md2.mesh.edges[:, 2] - 1] md2.mesh.edges[pos, 3] = Pelem[md2.mesh.edges[pos, 3] - 1] #remove edges when the 2 vertices are not in the domain. md2.mesh.edges = md2.mesh.edges[np.nonzero( np.logical_and(md2.mesh.edges[:, 0], md2.mesh.edges[:, 1]) )[0], :] #Replace all zeros by -1 in the last two columns pos = np.nonzero(md2.mesh.edges[:, 2] == 0)[0] md2.mesh.edges[pos, 2] = -1 pos = np.nonzero(md2.mesh.edges[:, 3] == 0)[0] md2.mesh.edges[pos, 3] = -1 #Invert -1 on the third column with last column (Also invert first two columns!!) pos = np.nonzero(md2.mesh.edges[:, 2] == -1)[0] md2.mesh.edges[pos, 2] = md2.mesh.edges[pos, 3] md2.mesh.edges[pos, 3] = -1 values = md2.mesh.edges[pos, 1] md2.mesh.edges[pos, 1] = md2.mesh.edges[pos, 0] md2.mesh.edges[pos, 0] = values #Finally remove edges that do not belong to any element pos = np.nonzero( np.logical_and(md2.mesh.edges[:, 1] == -1, md2.mesh.edges[:, 2] == -1))[0] md2.mesh.edges = np.delete(md2.mesh.edges, pos, axis=0) #Penalties if np.any(np.logical_not(np.isnan(md2.stressbalance.vertex_pairing))): for i in xrange(np.size(md1.stressbalance.vertex_pairing, axis=0)): md2.stressbalance.vertex_pairing[i, :] = Pnode[ md1.stressbalance.vertex_pairing[i, :]] md2.stressbalance.vertex_pairing = md2.stressbalance.vertex_pairing[ np.nonzero(md2.stressbalance.vertex_pairing[:, 0])[0], :] if np.any(np.logical_not(np.isnan(md2.masstransport.vertex_pairing))): for i in xrange(np.size(md1.masstransport.vertex_pairing, axis=0)): md2.masstransport.vertex_pairing[i, :] = Pnode[ md1.masstransport.vertex_pairing[i, :]] md2.masstransport.vertex_pairing = md2.masstransport.vertex_pairing[ np.nonzero(md2.masstransport.vertex_pairing[:, 0])[0], :] #recreate segments if md1.mesh.__class__.__name__ == 'mesh2d': md2.mesh.vertexconnectivity = NodeConnectivity( md2.mesh.elements, md2.mesh.numberofvertices)[0] md2.mesh.elementconnectivity = ElementConnectivity( md2.mesh.elements, md2.mesh.vertexconnectivity)[0] md2.mesh.segments = contourenvelope(md2) md2.mesh.vertexonboundary = np.zeros(numberofvertices2, bool) md2.mesh.vertexonboundary[md2.mesh.segments[:, 0:2] - 1] = True else: #First do the connectivity for the contourenvelope in 2d md2.mesh.vertexconnectivity = NodeConnectivity( md2.mesh.elements2d, md2.mesh.numberofvertices2d)[0] md2.mesh.elementconnectivity = ElementConnectivity( md2.mesh.elements2d, md2.mesh.vertexconnectivity)[0] segments = contourenvelope(md2) md2.mesh.vertexonboundary = np.zeros( numberofvertices2 / md2.mesh.numberoflayers, bool) md2.mesh.vertexonboundary[segments[:, 0:2] - 1] = True md2.mesh.vertexonboundary = np.tile(md2.mesh.vertexonboundary, md2.mesh.numberoflayers) #Then do it for 3d as usual md2.mesh.vertexconnectivity = NodeConnectivity( md2.mesh.elements, md2.mesh.numberofvertices)[0] md2.mesh.elementconnectivity = ElementConnectivity( md2.mesh.elements, md2.mesh.vertexconnectivity)[0] #Boundary conditions: Dirichlets on new boundary #Catch the elements that have not been extracted orphans_elem = np.nonzero(np.logical_not(flag_elem))[0] orphans_node = np.unique(md1.mesh.elements[orphans_elem, :]) - 1 #Figure out which node are on the boundary between md2 and md1 nodestoflag1 = np.intersect1d(orphans_node, pos_node) nodestoflag2 = Pnode[nodestoflag1].astype(int) - 1 if np.size(md1.stressbalance.spcvx) > 1 and np.size( md1.stressbalance.spcvy) > 2 and np.size( md1.stressbalance.spcvz) > 2: if np.size(md1.inversion.vx_obs) > 1 and np.size( md1.inversion.vy_obs) > 1: md2.stressbalance.spcvx[nodestoflag2] = md2.inversion.vx_obs[ nodestoflag2] md2.stressbalance.spcvy[nodestoflag2] = md2.inversion.vy_obs[ nodestoflag2] else: md2.stressbalance.spcvx[nodestoflag2] = np.nan md2.stressbalance.spcvy[nodestoflag2] = np.nan print "\n!! extract warning: spc values should be checked !!\n\n" #put 0 for vz md2.stressbalance.spcvz[nodestoflag2] = 0 if np.any(np.logical_not(np.isnan(md1.thermal.spctemperature))): md2.thermal.spctemperature[nodestoflag2] = 1 #Results fields if md1.results: md2.results = results() for solutionfield, field in md1.results.__dict__.iteritems(): if isinstance(field, list): setattr(md2.results, solutionfield, []) #get time step for i, fieldi in enumerate(field): if isinstance(fieldi, results) and fieldi: getattr(md2.results, solutionfield).append(results()) fieldr = getattr(md2.results, solutionfield)[i] #get subfields for solutionsubfield, subfield in fieldi.__dict__.iteritems( ): if np.size(subfield) == numberofvertices1: setattr(fieldr, solutionsubfield, subfield[pos_node]) elif np.size(subfield) == numberofelements1: setattr(fieldr, solutionsubfield, subfield[pos_elem]) else: setattr(fieldr, solutionsubfield, subfield) else: getattr(md2.results, solutionfield).append(None) elif isinstance(field, results): setattr(md2.results, solutionfield, results()) if isinstance(field, results) and field: fieldr = getattr(md2.results, solutionfield) #get subfields for solutionsubfield, subfield in field.__dict__.iteritems( ): if np.size(subfield) == numberofvertices1: setattr(fieldr, solutionsubfield, subfield[pos_node]) elif np.size(subfield) == numberofelements1: setattr(fieldr, solutionsubfield, subfield[pos_elem]) else: setattr(fieldr, solutionsubfield, subfield) #Keep track of pos_node and pos_elem md2.mesh.extractedvertices = pos_node + 1 md2.mesh.extractedelements = pos_elem + 1 return md2
def squaremesh(md, Lx, Ly, nx, ny): """ SQUAREMESH - create a structured square mesh This script will generate a structured square mesh Lx and Ly are the dimension of the domain (in meters) nx anx ny are the number of nodes in the x and y direction The coordinates x and y returned are in meters. Usage: [md]=squaremesh(md,Lx,Ly,nx,ny) """ #get number of elements and number of nodes nel = (nx - 1) * (ny - 1) * 2 nods = nx * ny #initialization index = np.zeros((nel, 3), int) x = np.zeros((nx * ny)) y = np.zeros((nx * ny)) #create coordinates for n in xrange(0, nx): for m in xrange(0, ny): x[n * ny + m] = float(n) y[n * ny + m] = float(m) #create index for n in xrange(0, nx - 1): for m in xrange(0, ny - 1): A = n * ny + (m + 1) B = A + 1 C = (n + 1) * ny + (m + 1) D = C + 1 index[n * (ny - 1) * 2 + 2 * m, :] = [A, C, B] index[n * (ny - 1) * 2 + 2 * (m + 1) - 1, :] = [B, C, D] #Scale x and y x = x / np.max(x) * Lx y = y / np.max(y) * Ly #create segments segments = np.zeros((2 * (nx - 1) + 2 * (ny - 1), 3), int) #left edge: segments[0:ny - 1, :] = np.vstack( (np.arange(2, ny + 1), np.arange(1, ny), (2 * np.arange(1, ny) - 1))).T #right edge: segments[ny - 1:2 * (ny - 1), :] = np.vstack( (np.arange(ny * (nx - 1) + 1, nx * ny), np.arange(ny * (nx - 1) + 2, nx * ny + 1), 2 * np.arange( (ny - 1) * (nx - 2) + 1, (nx - 1) * (ny - 1) + 1))).T #front edge: segments[2 * (ny - 1):2 * (ny - 1) + (nx - 1), :] = np.vstack( (np.arange(2 * ny, ny * nx + 1, ny), np.arange(ny, ny * (nx - 1) + 1, ny), np.arange(2 * (ny - 1), 2 * (nx - 1) * (ny - 1) + 1, 2 * (ny - 1)))).T #back edge segments[2 * (ny - 1) + (nx - 1):2 * (nx - 1) + 2 * (ny - 1), :] = np.vstack( (np.arange(1, (nx - 2) * ny + 2, ny), np.arange(ny + 1, ny * (nx - 1) + 2, ny), np.arange(1, 2 * (nx - 2) * (ny - 1) + 2, 2 * (ny - 1)))).T #plug coordinates and nodes md.mesh = mesh2d() md.mesh.x = x md.mesh.y = y md.mesh.numberofvertices = nods md.mesh.vertexonboundary = np.zeros((nods), bool) md.mesh.vertexonboundary[segments[:, 0:2] - 1] = True #plug elements md.mesh.elements = index md.mesh.segments = segments md.mesh.numberofelements = nel #Now, build the connectivity tables for this mesh. md.mesh.vertexconnectivity = NodeConnectivity(md.mesh.elements, md.mesh.numberofvertices)[0] md.mesh.elementconnectivity = ElementConnectivity( md.mesh.elements, md.mesh.vertexconnectivity)[0] return md