rbmatlab 0.10.01
Public Member Functions | Static Public Member Functions | Public Attributes
triagrid Class Reference

Detailed Description

A triangular conforming grid in two dimensions.

General geometrical information is stored including neighbour information. Therefore also boundary neighbour relations can be specified. Boundary edges can be marked by "painting rectangles": the edges with midpoints within such a rectangle are marked accordingly. By this boundary edges can be marked for later special use. Much (partially redundant) information is stored in the grid, which might be useful in simulations.

Definition at line 1 of file triagrid.m.

Inheritance diagram for triagrid:
Inheritance graph
[legend]
Collaboration diagram for triagrid:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 triagrid ( varargin)
 constructor of a triangular conform grid in 2 dimensions with following Synopsis:
function  demo ()
 small script demonstrating the possibilities of the triagrid class.
function  display ()
 display method for triagrid
function gridp gridpart ( eind)
 function extracting a part of a triagrid defined by the given element indices in the vector eind.
function lcoord llocal2local ( faceinds, llcoord)
 function performing a 1D edge-local coordinate (lcoord) to 2D local coordinate transformation of given faces
function glob local2global ( einds, loc, params)
 function performing a local to global coordinate change of vectors of coordinate pairs.
function p plot ( params)
 plot of a rectgrid via plot_polygon_grid()
function gridbase grid set_nbi ( nbind, values)
 function setting some neighbour indices of a grid to specified values.
function triagrid gcopy copy ()
 

returns a deep copy object of the grid implementation


Static Public Member Functions

static function [
C ,
G ] = 
aff_trafo_glob2loc ( x0, y0)
 function giving the coefficients for the affine transformation from original/global triangle to the reference/local one.
static function [
C ,
G ] = 
aff_trafo_loc2glob ( x0, y0)
 function giving the coefficients for the affine transformation from reference/local triangle to the original/global one.
static function [
C ,
G ] = 
aff_trafo_orig2ref ( x0, y0, varargin)
 function giving the coefficients for the affine transformation from original triangle to the reference one,
static function loc global2local ( gridbase grid, elementid, glob)
 function getting a triagrid, an element-ID and a vector of points and giving a vector of transformed points. The triangle given by elementid is used for the creation of an affine map to the standard tringle, then this transformation is used for all the points in glob
static function
micro2macro
micro2macro_map ( microgrid, macrogrid)
 function defining a vector micro2macro containing the information which triangle of the microgrid lies in which triangle of the macrogrid, defined in the model micro2macro(5) = 7 means that micro-triangle nr 5 lies in macro-triangle nr 7

Public Attributes

 nedges_interior
 number of interior edges
 nedges_boundary
 number of boundary edges

Constructor & Destructor Documentation

triagrid.triagrid (   varargin)

constructor of a triangular conform grid in 2 dimensions with following Synopsis:

  • triagrid() : construction of a default triagrid (2d unit square, loaded from file, -1 as outer neighbour indices)
  • triagrid(tgrid) : copy-constructor
  • triagrid(params) : in this case either
    1. the field grid_initfile is existent in params. Then the file is read and a pointlist p and a triangle. Procedure is then executing the following constructor
    2. a structured equidistant triangle grid is generated
      • params.xnumintervals : number of elements along x directions
      • params.ynumintervals : number of elements along y directions
      • params.xrange,yrange : intervals covered along the axes
      with the optional fields
      • params.bnd_rect_corner1 : coordinates of lower corner of to be marked boundaries
      • params.bnd_rect_corner2 : coordinates of upper corner of to be marked boundaries
      • params.bnd_rect_index: integer index to be set on the edges in the above defined rectangle. Should not be positive integer in the range of the number of elements. use negative indices for certain later discrimination.
      For the last three optional boundary settings, also multiple rectangles can be straightforwardly defined by accepting matrix of columnwise corners1, corners2 and a vector of indices for the different rectangles.
  • triagrid(p,t,params) : generate triagrid from triangle-data with certain options.
    • p is assumed to be a 2 x npoints matrix with coordinates
    • t is assumed to be a XX x ntriangles matrix, but only first three rows are used == vertex indices, in clockwise order as default, all nondefined edges are set to -1, then the following "rectangles" are set as specified in params
Using this class, grids from PDETOOLS can be used

pdetools => generate your grid and export p and t to MATLAB workspace

         save('mygrid','p','t')
         grid = triagrid(struct('grid_initfile','mygrid'));

perhaps later: constructor by duneDGF-file? perhaps later: contructor-flag: full vs non-full => only compute redundant information if required.

Note:
for diffusion-discretization with FV-schemes, points $S_i$ must exist, such that $S_i S_j$ is perpendicular to edge i j, the intersection points are denoted $S_{ij}$
Parameters:
vararginvariable number of input arguments, see above for description of possible configurations.
Generated fields of grid:
  • nelements —  number of elements
  • nvertices —  number of vertices
  • nneigh —  3
  • A —  vector of element area
  • Ainv —  vector of inverted element area
  • X —  vector of vertex x-coordinates
  • Y —  vector of vertex y-coordinates
  • VI —  matrix of vertex indices: VI(i,j) is the global index of j-th vertex of element i
  • CX —  vector of centroid x-values
  • CY —  vector of centroid y-values
  • NBI —  NBI(i,j) = element index of j-th neighbour of element i boundary faces are set to -1 or negative values are requested by params.boundary_type
  • INB —  INB(i,j) = local edge number in NBI(i,j) leading from element NBI(i,j) to element i, i.e. satisfying NBI(NBI(i,j), INB(i,j)) = i
  • EL —  EL(i,j) = length of edge from element i to neighbour j
  • DC —  DC(i,j) = distance from centroid of element i to NB j for boundary elements, this is the distance to the reflected element (for use in boundary treatment)
  • NX —  NX(i,j) = x-coordinate of unit outer normal of edge from el i to NB j
  • NY —  NY(i,j) = y-coordinate of unit outer normal of edge from el i to NB j
  • ECX —  ECX(i,j) = x-coordinate of midpoint of edge from el i to NB j
  • ECY —  ECY(i,j) = y-coordinate of midpoint of edge from el i to NB j
  • SX —  vector of x-coordinates of point $S_i$ (for rect: identical to centroids)
  • SY —  vector of y-coordinate of point $S_j$ (for rect: identical to centroids)
  • ESX —  ESX(i,j) = x-coordinate of point $S_{ij}$ on edge el i to NB j
  • ESY —  ESY(i,j) = y-coordinate of point $S_{ij}$ on edge el i to NB j
  • DS —  DS(i,j) = distance from $S_i$ of element i to $S_j$ of NB j for boundary elements, this is the distance to the reflected element (for use in boundary treatment)
  • hmin —  minimal element-diameter
  • alpha —  geometry bound (simultaneously satisfying $\alpha \cdot h_i^d \leq A(T_i)$, $\alpha \cdot \mbox{diam}(T_i) \leq h_i^{d-1}$ and $\alpha \cdot h_i \leq $dist(midpoint $i$ to any neigbour) )
  • JIT —  Jacobian inverse transposed JIT(i,:,:) is a 2x2-matrix of the Jac. Inv. Tr. on element i

Definition at line 39 of file triagrid.m.


Member Function Documentation

function [ C , G ] = triagrid.aff_trafo_glob2loc (   x0,
  y0 
) [static]

function giving the coefficients for the affine transformation from original/global triangle to the reference/local one.

In detail, this implements a transformation of the type

\[T_{i,aff}(x;\mu) = C_{i,aff}(\mu) + \sum_{j=1,2} G_{ij}^k(\mu) x_j \qquad i=1,2\]

     triangle:                     /| (x0(3),y0(3))              (0,1)  |\
                                  / |                ---T--->           | \
                   (x0(1),y0(1)) /__| (x0(2),y0(2))              (0,0)  |__\ (1,0)
     
See also:
aff_trafo_loc2glob() which gives the transformation in the other direction (local to global)
Parameters:
x0vector of size 3 x 1 holding x values of the original/global triangle
y0vector of size 3 x 1 holding y values of the original/global triangle
Return values:
Cmatrix of size 2 x 1 with entries C=[c1; c2]
Gmatrix of size 2 x 2 with entries G=[g11, g12; g21, g22]

Definition at line 2 of file aff_trafo_glob2loc.m.

function [ C , G ] = triagrid.aff_trafo_loc2glob (   x0,
  y0 
) [static]

function giving the coefficients for the affine transformation from reference/local triangle to the original/global one.

In detail, this implements a transformation of the type

\[T^{-1}_{i,aff}(x;\mu) = C^{-1}_{i,aff}(\mu) + \sum_{j=1,2} G^{k,-1}_{ij}(\mu) x_j \qquad i=1,2\]

The transformation is calculated for the standard triangle
     triangle:   (0,1)  |\                                   /| (x0(3),y0(3))
                        | \        ---T--->                 / |
                 (0,0)  |__\ (1,0)           (x0(1),y0(1)) /__| (x0(2),y0(2))
     
See also:
aff_trafo_glob2loc() which gives the transformation in the other direction (global to local)
Parameters:
x0vector of size 3 x 1 holding x values of the original/global triangle
y0vector of size 3 x 1 holding y values of the original/global triangle
Return values:
Cmatrix of size 2 x 1 with entries C=[c1; c2]
Gmatrix of size 2 x 2 with entries G=[g11, g12; g21, g22]

Definition at line 2 of file aff_trafo_loc2glob.m.

function [ C , G ] = triagrid.aff_trafo_orig2ref (   x0,
  y0,
  varargin 
) [static]

function giving the coefficients for the affine transformation from original triangle to the reference one,

\[T^{-1}_{i,aff}(x;\mu) = C^{-1}_{i,aff}(\mu) + \sum_{j=1,2} G^{k,-1}_{ij}(\mu) x_j \qquad i=1,2\]

     triangle:                     /| (x0(3),y0(3))              (0,1)  |\
                                  / |                ---T--->           | \
                   (x0(1),y0(1)) /__| (x0(2),y0(2))              (0,0)  |__\ (1,0)
     
Deprecated:
I guess, that this function is deprecated and aff_trafo_glob2loc() should be used instead...

function giving c1, c2, g11, g12, g21, g22 so: C=[c1; c2] and G=[g11, g12; g21, g22]

Parameters:
x0x0
y0y0
vararginvarargin
Return values:
CC
GG

Definition at line 2 of file aff_trafo_orig2ref.m.

function triagrid gcopy = triagrid.copy ( ) [virtual]

returns a deep copy object of the grid implementation

Return values:
gcopyobject This is a deep copy of the current instance.

Implements gridbase.

Definition at line 988 of file triagrid.m.

function triagrid.display ( )

display method for triagrid

This inherits gridbase.display() and adds information on

  • the number of interior edges and
  • the number of boundary edges.

Reimplemented from gridbase.

Definition at line 2 of file display.m.

function loc = triagrid.global2local ( gridbase  grid,
  elementid,
  glob 
) [static]

function getting a triagrid, an element-ID and a vector of points and giving a vector of transformed points. The triangle given by elementid is used for the creation of an affine map to the standard tringle, then this transformation is used for all the points in glob

input
grid: triagrid elementid: scalar nr of the triangle in triagrid, which defines the affine map glob: n-by-2-vector of points to be transformed
output
loc: n-by-2-vector of the transformed points

Oliver Zeeb, 02.02.11

Parameters:
gridan object
elementidelementid
globglob
Return values:
locloc
Required fields of grid:
  • VI —  matrix of vertex indices: VI(i,j) is the global index of the j-th vertex of element i
  • X —  vector of vertex $x$-coordinates.
  • Y —  vector of vertex $y$-coordinates.

Definition at line 2 of file global2local.m.

function gridbase gridp = triagrid.gridpart (   eind)

function extracting a part of a triagrid defined by the given element indices in the vector eind.

Note:
The neighbour information of the new resulting boundaries is set to -10
The properties gridbase.hmin, gridbase.alpha and the distance-information in the new boundary elements are simply copied. I.e. these fields do not completely meet the definition in the constructor. They might be chosen slightly different, such that the gridp would be really identical to the result generated by the constructor on the subset of points.
Parameters:
eindvector of cell indices which shall be extracted from the grid.
Return values:
gridpthe partial grid with extracted cells eind.
Generated fields of gridp:
  • nelements —  nelements
  • nvertices —  nvertices
  • A —  A
  • Ainv —  Ainv
  • X —  X
  • Y —  Y
  • VI —  VI
  • CX —  CX
  • CY —  CY
  • NBI —  NBI
  • INB —  INB
  • EL —  EL
  • DC —  DC
  • NX —  NX
  • NY —  NY
  • ECX —  ECX
  • ECY —  ECY
  • SX —  SX
  • SY —  SY
  • ESX —  ESX
  • ESY —  ESY
  • DS —  DS
  • JIT —  JIT

Reimplemented from gridbase.

Definition at line 2 of file gridpart.m.

function lcoord = triagrid.llocal2local (   faceinds,
  llcoord 
)

function performing a 1D edge-local coordinate (lcoord) to 2D local coordinate transformation of given faces

Parameters:
faceindsThe face indices on which the transformation shall take place. (1:3)
llcoordsingle real number between 0 and 1 defining the edge-local vertex coordinate.
Return values:
lcoordmatrix of size 2 x |faceinds| holding the local coordinates.

Definition at line 2 of file llocal2local.m.

function glob = triagrid.local2global (   einds,
  loc,
  params 
)

function performing a local to global coordinate change of vectors of coordinate pairs.

If the three vertices of a triangle are v1,v2,v3, then the global coordinate of a single point is

 glob = v1 + loc(:,1).*(v2-v1) + loc(:,2).*(v3-v1); 
Parameters:
eindsvector of cell indices $i_k$, $k=1,...,K$.
locmatrix of size $K \times 2$ holding local barycentric coordinate pairs for each cell index $i_k$, $k=1,...,K$.
paramsparams
Return values:
globglobal coordinate pairs [X, Y] with vectors X and Y of length $K$.

Definition at line 2 of file local2global.m.

function micro2macro = triagrid.micro2macro_map (   microgrid,
  macrogrid 
) [static]

function defining a vector micro2macro containing the information which triangle of the microgrid lies in which triangle of the macrogrid, defined in the model micro2macro(5) = 7 means that micro-triangle nr 5 lies in macro-triangle nr 7

microgrid and macrogrid must be triagrid

Oliver Zeeb, 01.02.11

Parameters:
microgridmicrogrid
macrogridmacrogrid
Return values:
micro2macromicro2macro
Required fields of microgrid:
  • nelements —  nelements
  • X —  X
  • Y —  Y
  • nvertices —  nvertices
  • VI —  VI
Required fields of macrogrid:
  • X —  X
  • Y —  Y
  • VI —  VI
  • nelements —  nelements

Definition at line 2 of file micro2macro_map.m.

function p = triagrid.plot (   params)

plot of a rectgrid via plot_polygon_grid()

see help plot_polygon_grid() for further information

plot method for a 2D polygonal grid. This routine can be used for triangular and rectangular grids.

A line plot is performed as default.

Todo:
For large grids, the routine can be slow. In these cases interestingly, the grid plotting should be implemented with patches, as that seems to be faster...
Parameters:
paramsoptional structure holding fields controlling the plot output.
Return values:
pThis is the list of handles to the graphics primitives
Required fields of params:
  • axis_tight —  axis tight
Optional fields of params:
  • color —  RGB vector of line/patch color
  • shrink_factor —  if this flag is given, the elements are plotted shrinked
  • plot_patch —  if this flag is set the plot is done by colored patches
  • axis_equal —  if this flag is set, set axis to equal scale
Parameters:
paramsparams
Return values:
pp

Definition at line 2 of file plot.m.

function gridbase grid = triagrid.set_nbi (   nbind,
  values 
)

function setting some neighbour indices of a grid to specified values.

Parameters:
nbindneighbour indices
valuesnew neighbour indices. This can be a single scalar or a vector of the same length as nbind.
Return values:
gridgrid

Definition at line 2 of file set_nbi.m.


The documentation for this class was generated from the following files:
All Classes Namespaces Files Functions Variables