![]() |
AdH Kraken 0.0.0
Next generation Adaptive Hydraulics
|
#include "adh.h"
Functions | |
void | get_cell_dofs_ivars (int *local_dofs, int **ivars, int nnodes, int *local_node_ids, int elem_nvars, int *elem_var_pos) |
Routine that gives an array of degrees of freedom local to the current process for a CG element using ivars double array Designed to work for nodal (CG) based dof mappings. More... | |
void | global_to_local_dbl_ivars (double *local_vals, int *nodeIDs, int nnodes, int *map_array, double *global_vals) |
This function takes a global array of values, and picks out the values for an array of nodes given a map array which takes nodeID -> dof index for the specific variable. (works for CG only) More... | |
void | global_to_local_SVECT2D_ivars (SVECT2D *local_vals, int *nodeIDs, int nnodes, int **ivars, int varx, int vary, double *global_vals) |
extracts sub-array of solution values for specific variable at a given set of nodes More... | |
This file collections functions responsible for finding order of finite element degrees of freedom for CG elements
void get_cell_dofs_ivars | ( | int * | local_dofs, |
int ** | ivars, | ||
int | nnodes, | ||
int * | local_node_ids, | ||
int | elem_nvars, | ||
int * | elem_var_pos | ||
) |
Routine that gives an array of degrees of freedom local to the current process for a CG element using ivars double array Designed to work for nodal (CG) based dof mappings.
[in,out] | local_dofs | (int*) - an array of integers that will give the degree of freedom numbers (equation numbers) for a given element local to the process |
[in] | fmaplocal | (int*) - an array of integers that gives the lowest d.o.f at a given node |
[in] | nnodes | (int) - the number of nodes on the element |
[in] | local_node_ids | (int*) - array of length nnodes containing node numbers local to process |
[in] | elem_nvars | (int) - number of solution variables active on the element |
[in] | elem_vars | (int*) - array of length elem_nvars that has the integer code for each variable |
[in] | node_physics_mat | (SMAT_PHYSICS*) - an array of SMAT_PHYSICS structs that contains variable info for each node |
[in] | nodal_physics_mat_id | (int*) - array of integers that gives the nodal physics mat id |
void global_to_local_dbl_ivars | ( | double * | local_vals, |
int * | nodeIDs, | ||
int | nnodes, | ||
int * | map_array, | ||
double * | global_vals | ||
) |
This function takes a global array of values, and picks out the values for an array of nodes given a map array which takes nodeID -> dof index for the specific variable. (works for CG only)
[in,out] | local_vals | (double*) - array of values corresponding to the node Ids |
[in] | nodeIDs | (int*) - node IDs (from the grid) |
[in] | nnodes | (int) - number of nodes requested to pull out from global array |
[in] | map_array | (int*) - array that takes node ID and produces dof # for a specific variable |
[in] | global_vals | (double*) - the global array of values that we want to pull from |
void global_to_local_SVECT2D_ivars | ( | SVECT2D * | local_vals, |
int * | nodeIDs, | ||
int | nnodes, | ||
int ** | ivars, | ||
int | varx, | ||
int | vary, | ||
double * | global_vals | ||
) |
extracts sub-array of solution values for specific variable at a given set of nodes
[in,out] | local | (double*) - the array containing local (to the element) values of a specific variable |
[in] | global | (double*) - the full array of the solution vector |
[in] | nodes | (int*) - array of node IDs local to process |
[in] | nnodes | (int) - the numer of nodes in the nodes array |
[in] | var | (int) - the variable code to be extracted |
[in] | fmaplocal | (int*) - array containing first dof number (local to process) at each node |
[in] | node_physics_mat | (SMAT_PHYSICS*) - array of SMAT_PHYSICS structs containing variable info at nodes |
[in] | nodal_physics_mat_id | (int*) - array of ints containing the nodal physics mat id at each node |