public static double[] projection(int fem_node_num, double[] fem_node_xy, int fem_element_order, int fem_element_num, int[] fem_element_node, int[] fem_element_neighbor, int fem_value_dim, double[] fem_value, int sample_node_num, double[] sample_node_xy) //****************************************************************************80 // // Purpose: // // PROJECTION evaluates an FEM function on a T3 or T6 triangulation. // // Discussion: // // Note that the sample values returned are true values of the underlying // finite element function. They are NOT produced by constructing some // other function that interpolates the data at the finite element nodes // (something which MATLAB's griddata function can easily do.) Instead, // each sampling node is located within one of the associated finite // element triangles, and the finite element function is developed and // evaluated there. // // MATLAB's scattered data interpolation is wonderful, but it cannot // be guaranteed to reproduce the finite element function corresponding // to nodal data. This routine can (or at least tries to//). // // So if you are using finite elements, then using THIS routine // (but not MATLAB's griddata function), what you see is what you have// // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 31 July 2009 // // Author: // // John Burkardt // // Parameters: // // Input, int FEM_NODE_NUM, the number of nodes. // // Input, double FEM_NODE_XY[2*FEM_NODE_NUM], the coordinates // of the nodes. // // Input, int FEM_ELEMENT_ORDER, the order of the elements, // either 3 or 6. // // Input, int FEM_ELEMENT_NUM, the number of triangles. // // Input, int FEM_ELEMENT_NODE(FEM_ELEMENT_ORDER,FEM_ELEMENT_NUM), the // nodes that make up each triangle. // // Input, int FEM_ELEMENT_NEIGHBOR[3*FEM_ELEMENT_NUM], the // index of the neighboring triangle on each side, or -1 if no neighbor there. // // Input, int FEM_VALUE_DIM, the "dimension" of the values. // // Input, double FEM_VALUE[FEM_VALUE_DIM*FEM_NODE_NUM], the // finite element coefficient values at each node. // // Input, int SAMPLE_NODE_NUM, the number of sample nodes. // // Input, double SAMPLE_NODE_XY[2*SAMPLE_NODE_NUM], the sample nodes. // // Output, double PROJECTION[FEM_VALUE_DIM*SAMPLE_NODE_NUM], // the sampled values. // { int edge = 0; int j; double[] p_xy = new double[2]; int t = 0; DelaunaySearchData data = new(); double alpha = 0; double beta = 0; double gammaa = 0; int step = 0; double[] b = new double[fem_element_order]; double[] dbdx = new double[fem_element_order]; double[] dbdy = new double[fem_element_order]; double[] sample_value = new double[fem_value_dim * sample_node_num]; int[] t_node = new int[fem_element_order]; double[] t_xy = new double[2 * fem_element_order]; // // For each sample point: find the triangle T that contains it, // and evaluate the finite element function there. // for (j = 0; j < sample_node_num; j++) { p_xy[0] = sample_node_xy[0 + 2 * j]; p_xy[1] = sample_node_xy[1 + 2 * j]; // // Find the triangle T that contains the point. // Search.triangulation_search_delaunay_a(ref data, fem_node_num, fem_node_xy, fem_element_order, fem_element_num, fem_element_node, fem_element_neighbor, p_xy, ref t, ref alpha, ref beta, ref gammaa, ref edge, ref step); switch (t) { case -1: Console.WriteLine(""); Console.WriteLine("PROJECTION - Fatal error!"); Console.WriteLine(" Triangulation search failed."); return(null); } // // Evaluate the finite element basis functions at the point in T. // int i; for (i = 0; i < fem_element_order; i++) { t_node[i] = fem_element_node[(i + t * fem_element_order) % fem_element_node.Length]; t_xy[0 + i * 2] = fem_node_xy[0 + t_node[i] * 2]; t_xy[1 + i * 2] = fem_node_xy[1 + t_node[i] * 2]; } switch (fem_element_order) { case 3: Basis_mn.basis_mn_t3(t_xy, 1, p_xy, ref b, ref dbdx, ref dbdy); break; case 6: Basis_mn.basis_mn_t6(t_xy, 1, p_xy, ref b, ref dbdx, ref dbdy); break; } // // Multiply by the finite element values to get the sample values. // for (i = 0; i < fem_value_dim; i++) { double dot = 0.0; int k; for (k = 0; k < fem_element_order; k++) { dot += fem_value[i + t_node[k] * fem_value_dim] * b[k]; } sample_value[i + j * fem_value_dim] = dot; } } return(sample_value); }
public static void derivative_average_t3(int node_num, double[] node_xy, int element_num, int[] element_node, double[] c, double[] dcdx, double[] dcdy) //****************************************************************************80 // // Purpose: // // DERIVATIVE_AVERAGE_T3 averages derivatives at the nodes of a T3 mesh. // // Discussion: // // This routine can be used to compute an averaged nodal value of any // quantity associated with the finite element function. At a node // that is shared by several elements, the fundamental function // U will be continuous, but its spatial derivatives, for instance, // will generally be discontinuous. This routine computes the // value of the spatial derivatives in each element, and averages // them, to make a reasonable assignment of a nodal value. // // Note that the ELEMENT_NODE array is assumed to be 1-based, rather // than 0-based. Thus, entries from this array must be decreased by // 1 before being used! // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 10 June 2006 // // Author: // // John Burkardt // // Parameters: // // Input, int NODE_NUM, the number of nodes. // // Input, double NODE_XY[2*NODE_NUM], the coordinates of the nodes. // // Input, int ELEMENT_NUM, the number of elements. // // Input, int ELEMENT_NODE[3*ELEMENT_NUM], the element->node data. // // Input, double C[NODE_NUM], the finite element coefficient vector. // // Output, double DCDX[NODE_NUM], DCDY[NODE_NUM], the averaged // values of dCdX and dCdY at the nodes. // { const int OFFSET = 1; double[] dphidx = new double[3 * 3]; double[] dphidy = new double[3 * 3]; int element; int node; int[] node_count = new int[node_num]; double[] phi = new double[3 * 3]; double[] t = new double[2 * 3]; for (node = 0; node < node_num; node++) { node_count[node] = 0; dcdx[node] = 0.0; dcdy[node] = 0.0; } // // Consider every element. // for (element = 0; element < element_num; element++) { // // Get the coordinates of the nodes of the element. // int j; for (j = 0; j < 3; j++) { int dim; for (dim = 0; dim < 2; dim++) { t[dim + 2 * j] = node_xy[dim + (element_node[j + element * 3] - OFFSET)]; } } // // Evaluate the X and Y derivatives of the 3 basis functions at the // 3 nodes. // Basis_mn.basis_mn_t3(t, 3, t, ref phi, ref dphidx, ref dphidy); // // Evaluate dCdX and dCdY at each node in the element, and add // them to the running totals. // int node_local1; for (node_local1 = 0; node_local1 < 3; node_local1++) { int node_global1 = element_node[node_local1 + element * 3] - OFFSET; int node_local2; for (node_local2 = 0; node_local2 < 3; node_local2++) { int node_global2 = element_node[node_local2 + element * 3] - OFFSET; dcdx[node_global1] += c[node_global2] * dphidx[node_local2 + node_local1 * 3]; dcdy[node_global1] += c[node_global2] * dphidy[node_local2 + node_local1 * 3]; } node_count[node_global1] += 1; } } // // Average the running totals. // for (node = 0; node < node_num; node++) { dcdx[node] /= node_count[node]; dcdy[node] /= node_count[node]; } }
public static double[] fem2d_evaluate(int fem_node_num, double[] fem_node_xy, int fem_element_order, int fem_element_num, int[] fem_element_node, int[] fem_element_neighbor, int fem_value_dim, double[] fem_value, int sample_node_num, double[] sample_node_xy) //****************************************************************************80 // // Purpose: // // GRID_SAMPLE samples a (scalar) FE function on a T3 or T6 triangulation. // // Discussion: // // Note that the values of U returned are true values of the underlying // finite element function. They are NOT produced by constructing some // other function that interpolates the data at the finite element nodes // (something which MATLAB's griddata function can easily do.) Instead, // each sampling node is located within one of the associated finite // element triangles, and the finite element function is developed and // evaluated there. // // MATLAB's scattered data interpolation is wonderful, but it cannot // be guaranteed to reproduce the finite element function corresponding // to nodal data. This routine can. // // So if you are using finite elements, then using THIS routine // (but not MATLAB's griddata function), what you see is what you have. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 01 June 2009 // // Author: // // John Burkardt // // Parameters: // // Input, int FEM_NODE_NUM, the number of nodes. // // Input, double FEM_NODE_XY[2*FEM_NODE_NUM], the coordinates of the nodes. // // Input, int FEM_ELEMENT_ORDER, the order of the triangles, either // 3 or 6. // // Input, int FEM_ELEMENT_NUM, the number of triangles. // // Input, int FEM_ELEMENT_NODE[FEM_ELEMENT_ORDER*FEM_ELEMENT_NUM], the // nodes that make up each element. // // Input, int FEM_ELEMENT_NEIGHBOR[3*FEM_ELEMENT_NUM], the index of the // neighboring triangle on each side, or -1 if no neighbor there. // // Input, int FEM_VALUE_DIM, the "dimension" of the values. // // Input, double FEM_VALUE[FEM_VALUE_DIM*FEM_NODE_NUM], the finite element // coefficient value at each node. // // Input, int SAMPLE_NODE_NUM, the number of sample nodes. // // Input, double SAMPLE_NODE_XY[2*SAMPLE_NODE_NUM], the sample nodes. // // Output, double SAMPLE_VALUE[FEM_VALUE_DIM*SAMPLE_NODE_NUM], // the sampled values. // { int edge = 0; int j; double[] p_xy = new double[2]; int t = 0; double[] b = new double[fem_element_order]; double[] dbdx = new double[fem_element_order]; double[] dbdy = new double[fem_element_order]; double[] sample_value = new double[fem_value_dim * sample_node_num]; int[] t_node = new int[fem_element_order]; double[] t_xy = new double[2 * fem_element_order]; // // For each sample point: find the triangle T that contains it, // and evaluate the finite element function there. // for (j = 0; j < sample_node_num; j++) { p_xy[0] = sample_node_xy[0 + j * 2]; p_xy[1] = sample_node_xy[1 + j * 2]; // // Find the triangle T that contains the point. // /* * Search.triangulation_search_delaunay(fem_node_num, fem_node_xy, * fem_element_order, fem_element_num, fem_element_node, * fem_element_neighbor, p_xy, ref t, ref edge); */ Search.triangulation_search_delaunay(fem_node_num, fem_node_xy, fem_element_order, fem_element_num, fem_element_node, fem_element_neighbor, p_xy, ref t, ref edge); // // Evaluate the finite element basis functions at the point in T. // // Note that in the following loop, we assume that FEM_ELEMENT_NODE // is 1-based, and that for convenience we can get away with making // T_NODE 0-based. // int order; for (order = 0; order < fem_element_order; order++) { t_node[order] = fem_element_node[order + (t - 1) * fem_element_order] - 1; } for (order = 0; order < fem_element_order; order++) { t_xy[0 + order * 2] = fem_node_xy[0 + t_node[order] * 2]; t_xy[1 + order * 2] = fem_node_xy[1 + t_node[order] * 2]; } switch (fem_element_order) { case 3: Basis_mn.basis_mn_t3(t_xy, 1, p_xy, ref b, ref dbdx, ref dbdy); break; case 6: Basis_mn.basis_mn_t6(t_xy, 1, p_xy, ref b, ref dbdx, ref dbdy); break; } // // Multiply by the finite element values to get the sample values. // int i; for (i = 0; i < fem_value_dim; i++) { double dot = 0.0; for (order = 0; order < fem_element_order; order++) { dot += fem_value[i + t_node[order]] * b[order]; } sample_value[i + j * fem_value_dim] = dot; } } return(sample_value); }