rbmatlab 0.10.01
|
00001 classdef gridbase < handle 00002 % Base class for all grid classes 00003 % 00004 00005 properties 00006 00007 nelements; % number of overall elements (leaf + nonleaf) 00008 00009 nvertices; % number of vertices 00010 00011 nneigh; % number of neighbor elements 00012 00013 A; % vector of element area 00014 00015 Ainv; % vector of inverted element area 00016 00017 % matrix of vertex indices: 'VI(i,j)' is the global index of 00018 % 'j'-th vertex of element 'i' 00019 VI; 00020 00021 X; % vector of vertex x-coordinates 00022 00023 Y; % vector of vertex y-coordinates 00024 00025 CX; % vector of centroid x-values 00026 00027 CY; % vector of centroid y-values 00028 00029 % 'NBI(i,j)' = element index of 'j'-th neighbour of element 'i' 00030 % boundary faces are set to -1 or negative values are requested 00031 % by 'params.boundary_type' 00032 NBI; 00033 00034 % 'INB(i,j)' = local edge number in 'NBI(i,j)' leading from 00035 % element 'NBI(i,j)' to element 'i', i.e. satisfying 00036 % 'NBI(NBI(i,j), INB(i,j)) = i' 00037 INB; 00038 00039 % 'EL(i,j)' = length of edge from element 'i' to neighbour 'j' 00040 EL; 00041 00042 % 'DC(i,j)' = distance from centroid of element 'i' to NB 'j' 00043 % for boundary elements, this is the distance to the reflected 00044 % element (for use in boundary treatment) 00045 DC; 00046 00047 % 'NX(i,j)' = x-coordinate of unit outer normal of edge from el 00048 % 'i' to NB 'j' 00049 NX; 00050 00051 % 'NY(i,j)' = y-coordinate of unit outer normal of edge from el 00052 % 'i' to NB 'j' 00053 NY; 00054 00055 % 'ECX(i,j)' = x-coordinate of midpoint of edge from el 'i' to 00056 % NB 'j' 00057 ECX; 00058 00059 % 'ECY(i,j)' = x-coordinate of midpoint of edge from el 'i' to 00060 % NB 'j' 00061 ECY; 00062 00063 % vector of x-coordinates of point `S_i` (for rect: identical to 00064 % centroids) 00065 % 00066 % for diffusion-discretization with FV-schemes, points `S_i` 00067 % must exist, such that `S_i S_j` is perpendicular to edge `i j` 00068 % the intersection are denoted `S_ij` 00069 SX; 00070 00071 % vector of y-coordinates of point `S_i` (for rect: identical to 00072 % centroids) 00073 % 00074 % for diffusion-discretization with FV-schemes, points `S_i` 00075 % must exist, such that `S_i S_j` is perpendicular to edge `i j` 00076 % the intersection are denoted `S_ij` 00077 SY; 00078 00079 % 'ESX(i,j)' = x-coordinate of point `S_ij` on edge el 'i' to 00080 % NB 'j' 00081 ESX; 00082 00083 % 'ESY(i,j)' = x-coordinate of point `S_ij` on edge el 'i' to 00084 % NB 'j' 00085 ESY; 00086 00087 % 'DS(i,j)' = distance from `S_i` of element 'i' to `S_j` of NB 00088 % 'j' for boundary elements, this is the distance to the 00089 % reflected element (for use in boundary treatment) 00090 DS; 00091 00092 hmin; % minimal element-diameter 00093 00094 % geometry bound (simultaneously satisfying 00095 % ``\alpha h_i^d \leq A(T_i),`` 00096 % ``\alpha \text{diam}(T_i) \leq h_i^{d-1}`` and 00097 % ``\alpha h_i \leq \text{distance(midpoint i to any neighbour)}`` 00098 alpha; 00099 00100 % Jacobian inverse transposed 'JIT(i,:,:)' is a '2x2'-matrix of 00101 % the Jacobian Inverse Transposed on element 'i' 00102 JIT; 00103 00104 end 00105 00106 methods 00107 00108 function obj = gridbase 00109 % default constructor for gridbase 00110 end 00111 00112 00113 res = check_consistency(grid); 00114 00115 display(grid); 00116 00117 F = edge_quad_eval(grid,elids,edgeids,degree,FF); 00118 00119 F = edge_quad_eval_mean(grid,elids,edgeids,degree,FF); 00120 00121 PP = edge_quad_points(grid,elids,edgeids,degree); 00122 00123 [P1, P2] = get_edge_points(grid,elids,edgeids); 00124 00125 ENBI = get_enbi(grid, edge, tstep); 00126 00127 gridp=gridpart(grid, eind); 00128 00129 inspect(grid, params); 00130 00131 p = plot_polygon_grid(grid,params); 00132 00133 p = plot_element_data(grid,data,plot_params); 00134 00135 plot_element_data_sequence(grid,data,plot_params); 00136 00137 p = plot_vertex_data(grid,data,params); 00138 00139 grid = set_boundary_types(grid,params); 00140 00141 end 00142 00143 methods(Abstract) 00144 00145 gcopy = copy(grid); 00146 % returns a deep copy object of the grid implementation 00147 % 00148 % Every grid implementations needs to provide a deep copy method, because 00149 % it is derived from a handle class. 00150 % 00151 % Return values: 00152 % gcopy: Deep copy of type gridbase. 00153 00154 end 00155 00156 end