NODAL_VARS_STRUCT * copy_and_rearrange_nv(NODAL_VARS_STRUCT *nv_old) /***************************************************************** * * copy_and_rearrange_nv(): * * This function copies an old nodal variable structure into * a new one, rearranging the order of the variables as it * goes. * This function sets the order in the solution vector. Right * now, the order of the solution vector is set to its original * GOMA one: * Outer to Inner: * var_type * sub_var_type * matID * This is not necessarily how it should be for optimal * speed, however. In particular the species unknowns should * be grouped together. Therefore, rearrangement of the * solution vector might occur in the future. *****************************************************************/ { int var_type, subvarIndex = 0; int *rec_used; NODAL_VARS_STRUCT *nv_new = alloc_struct_1(NODAL_VARS_STRUCT, 1); nv_new->list_index = nv_old->list_index; /* * Note: It is permissible to have nodes with zero unknowns defined * on them. For example, if the element is quadratic and the * variable interpolations are all Q1, the extra nodes in each * element will not have any unknowns assigned to them. This * is a waste, but it is permissible. */ if (nv_old->Num_Var_Desc == 0) { return nv_new; } rec_used = alloc_int_1(nv_old->Num_Var_Desc, 0); for (var_type = V_FIRST; var_type < V_LAST; var_type++) { if (var_type == MASS_FRACTION) { for (subvarIndex = 0; subvarIndex < upd->Max_Num_Species_Eqn; subvarIndex++) { rearrange_mat_increasing(var_type, subvarIndex, nv_old, nv_new, rec_used); } } else { rearrange_mat_increasing(var_type, 0, nv_old, nv_new, rec_used); } } safer_free((void **) &rec_used); return nv_new; }
void ReduceBcast_BOR(int *ivec, int length) /******************************************************************** * * ReduceBcast_BOR() * * Does a logical bitwize OR operation on a vector of ints * distibuted across different processors. ********************************************************************/ { #ifdef PARALLEL int k, err, *ivec_recv; if (length <= 0) { printf(" ReduceBcast_BOR Warning, length = %d\n", length); return; } if (ivec == NULL) { EH(-1, " ReduceBcast_BOR fatal Interface error"); } ivec_recv = alloc_int_1(length, INT_NOINIT); err = MPI_Reduce(ivec, ivec_recv, length, MPI_INT, MPI_BOR, 0, MPI_COMM_WORLD); if (err != MPI_SUCCESS) { EH(-1, " ReduceBcast_BOR fatal MPI Error"); } err = MPI_Bcast(ivec_recv, V_LAST, MPI_INT, 0, MPI_COMM_WORLD); if (err != MPI_SUCCESS) { EH(-1, " ReduceBcast_BOR fatal MPI Error"); } for (k = 0; k < V_LAST; k++) { ivec[k] = ivec_recv[k]; } safer_free((void **) &ivec_recv); #endif }
void build_elem_elem(Exo_DB *exo) { int ce; int count; int e; int ebi; int elem; int ename; int face; //int his_dim, her_dim; int i, j; int index; int len_curr; int len_prev; int len_intr; int length, length_new, num_faces; int iel; int ioffset; int n; int neighbor_name = -1; int node; int num_elem_sides; int num_nodes; int snl[MAX_NODES_PER_SIDE]; /* Side Node List - NOT Saturday Night Live! */ char err_msg[MAX_CHAR_ERR_MSG]; /* * Integer arrays used to find intersection sets of node->element lists. */ int prev_set[MAX_EPN]; /* list of elements attached to previous node*/ int curr_set[MAX_EPN]; /* list of elements attached to "this" node */ int interset[MAX_EPN]; /* values of hits between */ int ip[MAX_EPN]; /* indeces of hits for prev_set[] */ int ic[MAX_EPN]; /* indeces of hits for curr_set[] */ /* * If the element->node and node->element connectivities have not been * built, then we won't be able to do this task. */ if ( ! exo->elem_node_conn_exists || ! exo->node_elem_conn_exists ) { EH(-1, "Build elem->node before node->elem."); return; } /* * The number of elements connected via conventional faces may be deduced * from the number of elements and their type. */ exo->elem_elem_pntr = (int *) smalloc((exo->num_elems+1)*sizeof(int)); length = 0; for ( i=0; i<exo->num_elem_blocks; i++) { length += exo->eb_num_elems[i] * get_num_faces(exo->eb_elem_type[i]); } exo->elem_elem_list = (int *) smalloc(length*sizeof(int)); /* * Initialize... */ for ( i=0; i<length; i++) { exo->elem_elem_list[i] = UNASSIGNED_YET; } /* elem = 0; for ( ebi=0; ebi<exo->num_elem_blocks; ebi++) { num_elem_sides = get_num_faces(exo->eb_elem_type[ebi]); for ( e=0; e<exo->eb_num_elems[ebi]; e++) { exo->elem_elem_pntr[elem] = count; elem++; count += num_elem_sides; } } */ /* * Walk through the elements, block by block. */ count = 0; elem = 0; for ( ebi=0; ebi<exo->num_elem_blocks; ebi++) { num_elem_sides = get_num_faces(exo->eb_elem_type[ebi]); for ( e=0; e<exo->eb_num_elems[ebi]; e++,elem++) { exo->elem_elem_pntr[elem] = count; count += num_elem_sides; /* * Look at each side of the element, collecting a unique * list of integers corresponding to the minimum number of nodes * needed to identify an entire side. * * Typically, the same number of nodes as space dimensions are * needed, with exceptions being the various "sides" of shells, * beams and trusses... */ for ( face=0; face<num_elem_sides; face++) { /* * Given the element and the face construct the * list of node numbers that determine that side. */ /* * Later, we might not need *all* the nodes on a side, * particularly for high order elements. It may suffice * to check only as many nodes as space dimensions that * the element lives in... */ num_nodes = build_side_node_list(elem, face, ebi, exo, snl); #ifdef DEBUG fprintf(stderr, "Elem %d, face %d has %d nodes: ", elem, face, num_nodes); for ( i=0; i<num_nodes; i++) { fprintf(stderr, " %d", snl[i]); } fprintf(stderr, "\n"); #endif /* DEBUG */ /* * Cross check: for each node in the side there is a list * of elements connected to it. Beginning with all the * elements connected to the first node (except for this given * element), cross check with all the elements connected with * the 2nd node to build an intersection set of one element. */ for ( i=0; i<MAX_EPN; i++) { prev_set[i] = -1; curr_set[i] = -1; interset[i] = -1; } len_prev = 0; len_curr = 0; len_intr = 0; for ( n=0; n<num_nodes; n++) { /* * Copy this node's element list into a clean "curr_set" array * that will be intersected with any previously gathered * lists of elements that qualify as promiscuously in * contact with nodes... */ node = snl[n]; for ( i=0; i<MAX_EPN; i++) { curr_set[i] = -1; } len_curr = 0; #ifdef DEBUG fprintf(stderr, "Traversing n->e connectivity of node %d\n", node); #endif /* DEBUG */ for ( ce=exo->node_elem_pntr[node]; ce<exo->node_elem_pntr[node+1]; ce++) { ename = exo->node_elem_list[ce]; #ifdef DEBUG fprintf(stderr, "\telem %d\n", ename); #endif /* DEBUG */ /* * Go ahead and accumulate the self element name * just as a consistency check.... */ /* if ( ename != e ) { } */ /* * PKN: The current Goma use of ->elem_elem... * is such that this connectivity should list * connections like QUAD-BAR or HEX-SHELL. * So, I'll add this dimension matching conditional */ /* PRS (Summer 2012): Need to change this for shell stacks which have the same dim*/ /* We need however to consider a special case (as of 8/30/2012 * this is a SHELL-on-SHELL stack. Viz. two materials, each a shell material * which share not a side but a face. Since faces of shells are sides * in patran speak, we need some special logic. We need to avoid adding * the friend shell element (neighboring material) to the current shell element * even though each material has the same number of sides. * Here goes (BTW, I cannot find max-nodes-per-element anywhere!!!!) */ int shell_on_shell = 0; int flippy_flop = 0; int nbr_ebid; int nbr_num_elem_sides; nbr_ebid = fence_post(ename, exo->eb_ptr, exo->num_elem_blocks+1); EH(nbr_ebid, "Bad element block ID!"); nbr_num_elem_sides = get_num_faces(exo->eb_elem_type[nbr_ebid]); shell_on_shell = 0; flippy_flop = 0; if (exo->eb_id[ebi] < 100 && exo->eb_id[nbr_ebid] >= 100) flippy_flop=1; if (exo->eb_id[ebi] >= 100 && exo->eb_id[nbr_ebid] < 100) flippy_flop=1; if ((nbr_ebid != ebi) && (strstr(exo->eb_elem_type[nbr_ebid], "SHELL")) && (strstr(exo->eb_elem_type[ebi], "SHELL")) && flippy_flop) shell_on_shell = 1; // his_dim = elem_info(NDIM, exo->eb_elem_itype[ebi]); // her_dim = elem_info(NDIM, exo->eb_elem_itype[exo->elem_eb[ename]]); // if( his_dim == her_dim ) if (nbr_num_elem_sides == num_elem_sides && !shell_on_shell) { curr_set[len_curr] = ename; len_curr++; } } /* * The first node is special - we'll just compare * it with itself by making the "previous set" just the * same as the current set... */ if ( n == 0 ) { for ( i=0; i<MAX_EPN; i++) { prev_set[i] = curr_set[i]; } len_prev = len_curr; } #ifdef DEBUG fprintf(stderr, "\ncurr_set: "); for ( i=0; i<len_curr; i++) { fprintf(stderr, "%d ", curr_set[i]); } fprintf(stderr, "\nprev_set: "); for ( i=0; i<len_prev; i++) { fprintf(stderr, "%d ", prev_set[i]); } #endif /* DEBUG */ /* * First, clean the intersection list and the list of * hit indeces in the previous and current lists. * * Then find the intersection of the previous and current * sets of elements attached to the previous and current * nodes... */ for ( i=0; i<MAX_EPN; i++) { interset[i] = -1; ip[i] = -1; ic[i] = -1; } len_intr = 0; len_intr = int_intersect(prev_set, curr_set, len_prev, len_curr, ip, ic); #ifdef DEBUG fprintf(stderr, "num_hits = %d\n", len_intr); #endif /* DEBUG */ /* * Now, let's make the intersection set the next previous * set of elements, a standard for comparison. We should * eventually boil down to either one or zero elements * that qualify... */ for ( i=0; i<MAX_EPN; i++) { prev_set[i] = -1; } for ( i=0; i<len_intr; i++) { prev_set[i] = curr_set[ic[i]]; } len_prev = len_intr; } #ifdef DEBUG fprintf(stderr, "Element [%d], face [%d], local_node [%d]\n", elem, face, n); fprintf(stderr, "Intersection set length = %d\n", len_intr); #endif /* DEBUG */ /* * Now consider the different cases. */ if ( len_intr == 2 ) { /* * The boiled list contains self and one other element. */ if ( prev_set[0] == elem ) { neighbor_name = prev_set[1]; } else { neighbor_name = prev_set[0]; if ( prev_set[1] != elem ) { sr = sprintf(err_msg, "2 elems ( %d %d ) 1 should be %d!", prev_set[0], prev_set[1], elem); EH(-1, err_msg); } } } else if ( len_intr == 1 && prev_set[0] == elem ) { /* * The boiled list has one member, this element. * * The face must connect either to outer space or to * another processor. */ if ( Num_Proc == 1 ) { neighbor_name = -1; } else { neighbor_name = -1; /* * I am going to punt for now. Later, revisit this * condition and insert code to check for neighbor * processors containing all the same face nodes. * * EH(-1, "Not done yet..."); * */ /* * Check if ALL the nodes on this face belong * to another processors list of nodes. I.e., the * node must all be in the external node list of * and belong to the same external processor. */ } } /* * Pathological cases that normally should not occur.... */ else if ( len_intr == 0 ) { sr = sprintf(err_msg, "Elem %d, face %d should self contain!", elem, face); EH(-1, err_msg); } else if ( len_intr == 1 && prev_set[0] != elem ) { sr = sprintf(err_msg, "Elem %d, face %d only connects with elem %d ?", elem, face, prev_set[0]); EH(-1, err_msg); } else { sr = sprintf(err_msg, "Unknown elem-elem connection elem %d, face %d, len_intr=%d", elem, face, len_intr); WH(-1, err_msg); } /* * Now we know how to assign the neighbor name for this face * of the element. */ index = exo->elem_elem_pntr[elem] + face; exo->elem_elem_list[index] = neighbor_name; } /* end face loop this elem */ } /* end elem loop this elemblock */ } /* end elem block loop */ exo->elem_elem_pntr[exo->num_elems] = count; /* last fencepost */ exo->elem_elem_conn_exists = TRUE; if (Linear_Solver == FRONT) { /* * Now that we have elem_elem_pntr and elem_elem_list for our parallel * world, we are going to use them also for optimal element bandwidth * reduction ordering. We will use METIS, but METIS requires the CSR * format, which is compressed, viz. we need to remove the -1s. Here * we go */ /* First check for the assumption that all blocks have same number of element faces. Stop if they don't and issue an error to the next aspiring developer */ for ( i=0; i<exo->num_elem_blocks; i++) { if(get_num_faces(exo->eb_elem_type[0]) != get_num_faces(exo->eb_elem_type[i]) ) { EH(-1,"Stop! We cannot reorder these elements with METIS with elemement type changes"); } } /* Now begin */ exo->elem_elem_xadj = (int *) smalloc((exo->num_elems+1)*sizeof(int)); /*initialize */ for(e=0; e<exo->num_elems+1 ; e++) { exo->elem_elem_xadj[e] = exo->elem_elem_pntr[e]; } /* Recompute length of adjacency list by removing external edges */ length_new = 0; for (i = 0; i < length; i++) { if(exo->elem_elem_list[i] != -1) length_new++; } exo->elem_elem_adjncy = alloc_int_1(length_new, -1); /* Now convert */ ioffset=0; for(iel = 0; iel < exo->num_elems; iel++) { /* Big assumption here that all blocks have the same */ /* element type. Can be furbished later since this is */ /* just for the frontal solver */ num_faces = get_num_faces(exo->eb_elem_type[0]); for (i= iel*num_faces; i < (iel+1)*num_faces; i++) { j = i - ioffset; if(exo->elem_elem_list[i] == -1) { ioffset++; for(e=iel +1; e <exo->num_elems+1; e++)exo->elem_elem_xadj[e]--; } else { exo->elem_elem_adjncy[j] = exo->elem_elem_list[i]; } } } /* convert to Fortran style */ for(e=0; e<exo->num_elems+1 ; e++) exo->elem_elem_xadj[e]++; for ( i=0; i<length_new; i++) exo->elem_elem_adjncy[i]++; } /* End FRONTAL_SOLVER if */ /* * Verification that every element/face has assigned something besides * the initial default value of "unassigned". * * For your convenience - FORTRAN 1-based numbering. */ #ifdef DEBUG for ( e=0; e<exo->num_elems; e++) { fprintf(stdout, "Elem %3d:", e+1); for ( ce=exo->elem_elem_pntr[e]; ce<exo->elem_elem_pntr[e+1]; ce++) { if ( exo->elem_elem_list[ce] == -1 ) { fprintf(stdout, " spc"); } else if ( exo->elem_elem_list[ce] < -1 ) { fprintf(stdout, " prc"); } else { fprintf(stdout, " %3d", exo->elem_elem_list[ce] + 1); } if ( exo->elem_elem_list[ce] == UNASSIGNED_YET ) { sr = sprintf(err_msg, "You need to plug a leak at elem (%d) face (%d)", exo->elem_elem_list[ce] + 1, ce - exo->elem_elem_pntr[e] + 1); EH(-1, err_msg); } } fprintf(stdout, "\n"); } #endif /* DEBUG */ #if FALSE demo_elem_elem_conn(exo); #endif return; }
void bc_matrl_index(Exo_DB *exo) /******************************************************************** * * bc_matrl_index(): * * Find out what materials are on each side of a boundary * condition. Note, some boundary conditions require one * to specify this information in the form of an element block * id. However, some others do not. This procedure attempts to * calculate this for all boundary conditions and then fill in * the BC_matrl_index_# elements of the Boundary_Condition * structure. * * The basic algorithm involves finding some information about * the side or node set on which the boundary condition is * applied. Then, given the bc name and this information, * a decision is made as to what the BC_matrl_index_#'s should be * set to. * * We gather the following information to make this decision: * 1) element block indecises input from deck * 2) Number of nodes in the bc set in each material * 3) Number of elements, whoses sides are in the side set, * in each material * 4) Node in the bc set which contains the minimum * number of materials * 5) Node in the bc set which contains the maximum * number of materials * 6) A random node in the bc set which contains a * specific number of materials (1, 2, 3, 4) * * * HKM NOTE: * This routine is a work in progress. Calculation of EDGE's and * VERTICES are not done yet. Also, mp aspects haven't been * figured out yet. * *******************************************************************/ { int ibc, ss_index, side_index, k, node_num, i; int i_apply_meth, num_matrl_needed = -1; int found = FALSE, min_node_matrl, max_node_matrl, min_matrl, max_matrl, node_matrl_1, node_matrl_2, node_matrl_3, node_matrl_4, matrl_first; int *bin_matrl, *ind_matrl, *bin_matrl_elem, *node_flag_1ss, node_count; int mat_index, success, ielem; NODE_INFO_STRUCT *node_ptr; UMI_LIST_STRUCT *matrlLP; struct Boundary_Condition *bc_ptr; static char *yo = "bc_matrl_index :"; bin_matrl = alloc_int_1(upd->Num_Mat * 4, INT_NOINIT); ind_matrl = bin_matrl + upd->Num_Mat; bin_matrl_elem = ind_matrl + upd->Num_Mat; node_flag_1ss = alloc_int_1(exo->num_nodes, INT_NOINIT); for (ibc = 0; ibc < Num_BC; ibc++) { bc_ptr = BC_Types + ibc; found = FALSE; /* * Some boundary condition specifications already require * you to specify the element blocks on either side of the * side set. */ switch (bc_ptr->BC_Name) { /* * For these boundary conditions, the element block ID numbers * are in the firest two integer slots */ case POROUS_PRESSURE_BC: case DARCY_CONTINUOUS_BC: case Y_DISCONTINUOUS_BC: case POROUS_GAS_BC: case VP_EQUIL_BC: case VN_POROUS_BC: case FLUID_SOLID_BC: case SOLID_FLUID_BC: case NO_SLIP_BC: case FLUID_SOLID_CONTACT_BC: case SOLID_FLUID_CONTACT_BC: case T_CONTACT_RESIS_BC: case T_CONTACT_RESIS_2_BC: case LIGHTP_JUMP_BC: case LIGHTM_JUMP_BC: case LIGHTP_JUMP_2_BC: case LIGHTM_JUMP_2_BC: bc_ptr->BC_matrl_index_1 = map_mat_index(bc_ptr->BC_Data_Int[0]); bc_ptr->BC_matrl_index_2 = map_mat_index(bc_ptr->BC_Data_Int[1]); break; /* * For this boundary condition the element block numbers * are in the second and third integer slots */ case VL_EQUIL_BC: case YFLUX_DISC_RXN_BC: case DISCONTINUOUS_VELO_BC: bc_ptr->BC_matrl_index_1 = map_mat_index(bc_ptr->BC_Data_Int[1]); bc_ptr->BC_matrl_index_2 = map_mat_index(bc_ptr->BC_Data_Int[2]); break; } /* * Initialize quantities */ for (i = 0; i < 4 * upd->Num_Mat; i++) bin_matrl[i] = 0; for (i = 0; i < exo->num_nodes; i++) node_flag_1ss[i] = 0; max_node_matrl = min_node_matrl = -1; min_matrl = 4000000; max_matrl = -1; node_matrl_1 = node_matrl_2 = node_matrl_3 = node_matrl_4 = -1; /* * Determine how many materials each boundary condition needs * based on the type of the boundary condition */ i_apply_meth = BC_Types[ibc].desc->i_apply; if (i_apply_meth == CROSS_PHASE_DISCONTINUOUS || i_apply_meth == CROSS_PHASE) { num_matrl_needed = 2; } else if (i_apply_meth == SINGLE_PHASE) { num_matrl_needed = 1; } /* * Loop over the nodes in the side set looking up what * materials are located at each node * -> this will work for side sets. Need to do node sets * as well. */ if (!strcmp(bc_ptr->Set_Type, "SS")) { for (ss_index = 0; ss_index < exo->num_side_sets; ss_index++) { /* * This logic works for one side set specifications. For * two side set specifications (EDGES), we will have to go * with a calculation of the union of side sets. */ if (bc_ptr->BC_ID == exo->ss_id[ss_index]) { found = TRUE; for (side_index = 0; side_index < exo->ss_num_sides[ss_index]; side_index++) { /* * Locate the element number, find the material index, * then bin the result. */ ielem = exo->ss_elem_list[exo->ss_elem_index[ss_index]+side_index]; mat_index = find_mat_number(ielem, exo); bin_matrl_elem[mat_index]++; for (k = exo->ss_node_side_index[ss_index][side_index]; k < exo->ss_node_side_index[ss_index][side_index+1]; k++) { node_num = exo->ss_node_list[ss_index][k]; if (!node_flag_1ss[node_num]) { node_flag_1ss[node_num] = 1; node_ptr = Nodes[node_num]; matrlLP = &(node_ptr->Mat_List); /* * Bin the materials at this node for later usage. */ for (i = 0; i < matrlLP->Length; i++) { #ifdef DEBUG_IGNORE_ELEMENT_BLOCK_CAPABILITY if (matrlLP->List[i] < 0) { fprintf(stderr,"Material list contains negative number\n"); EH(-1,"logic error in ignoring an element block"); } #endif bin_matrl[matrlLP->List[i]]++; } /* * Find the max and min number of materials for a * node in this side set */ if (matrlLP->Length > max_matrl) { max_matrl = matrlLP->Length; max_node_matrl = node_num; } if (matrlLP->Length < min_matrl) { min_matrl = matrlLP->Length; min_node_matrl = node_num; } /* * Find representative nodes with specific * numbers of materials */ if (matrlLP->Length == 1) node_matrl_1 = node_num; if (matrlLP->Length == 2) node_matrl_2 = node_num; if (matrlLP->Length == 3) node_matrl_3 = node_num; if (matrlLP->Length == 4) node_matrl_4 = node_num; } } } } } /* End of side set loop */ } /* if SS */ /* * Node Sets */ if (!strcmp(bc_ptr->Set_Type, "NS")) { for (ss_index = 0; ss_index < exo->num_node_sets; ss_index++) { /* * This logic works for one side set specifications. For * two side set specifications (EDGES), we will have to go * with a calculation of the union of side sets. */ if (bc_ptr->BC_ID == exo->ns_id[ss_index]) { found = TRUE; /* * Loop over the number of nodes */ for (k = 0; k < exo->ns_num_nodes[ss_index]; k++) { node_num = exo->ns_node_list[exo->ns_node_index[ss_index]+k]; if (!node_flag_1ss[node_num]) { node_flag_1ss[node_num] = 1; node_ptr = Nodes[node_num]; matrlLP = &(node_ptr->Mat_List); /* * Bin the materials at this node for later usage. */ for (i = 0; i < matrlLP->Length; i++) { bin_matrl[matrlLP->List[i]]++; } /* * Find the max and min number of materials for a * node in this side set */ if (matrlLP->Length > max_matrl) { max_matrl = matrlLP->Length; max_node_matrl = node_num; } if (matrlLP->Length < min_matrl) { min_matrl = matrlLP->Length; min_node_matrl = node_num; } /* * Find representative nodes with specific * numbers of materials */ if (matrlLP->Length == 1) node_matrl_1 = node_num; if (matrlLP->Length == 2) node_matrl_2 = node_num; if (matrlLP->Length == 3) node_matrl_3 = node_num; if (matrlLP->Length == 4) node_matrl_4 = node_num; } } } } /* End of side set loop */ } /* if NS */ /* * Ok, we have obtained statistics on the side and node sets * let's make a decision */ if (found) { matrl_first = find_next_max(bin_matrl, ind_matrl, upd->Num_Mat); success = assign_matrl_2(bc_ptr, matrl_first); if (success < 0) { printf("%s P_%d: problem in assigning first matrl index in ibc %d, %d:\n", yo, ProcID, ibc, matrl_first); bc_matrl_index_print(bc_ptr, bin_matrl, bin_matrl_elem, min_node_matrl, max_node_matrl, node_matrl_1, node_matrl_2, node_matrl_3, node_matrl_4, ibc); } matrl_first = find_next_max(bin_matrl, ind_matrl, upd->Num_Mat); if (matrl_first >= 0) { success = assign_matrl_2(bc_ptr, matrl_first); if (success < 0) { printf("%s P_%d: problem in assigning second matrl index in ibc %d, %d:\n", yo, ProcID, ibc, matrl_first); bc_matrl_index_print(bc_ptr, bin_matrl, bin_matrl_elem, min_node_matrl, max_node_matrl, node_matrl_1, node_matrl_2, node_matrl_3, node_matrl_4, ibc); } } else { if (num_matrl_needed > 1) { printf("%s P_%d: problem in finding a needed second matrl index:\n", yo, ProcID); EH(-1,"bc_matrl_index ERROR"); } } } /* * For debug purposes, print out everything that we have found * out and decided about this bc on all of the processors. */ #ifdef DEBUG_HKM print_sync_start(FALSE); bc_matrl_index_print(bc_ptr, bin_matrl, bin_matrl_elem, min_node_matrl, max_node_matrl, node_matrl_1, node_matrl_2, node_matrl_3, node_matrl_4, ibc); print_sync_end(FALSE); #endif /* * MP Fix: We may not get the same results on all processors * In this case, just take the processor with the most * nodes in this bc and with a valid result, and use * that. Broadcast that result to all nodes. Cross your * fingers. */ node_count = 0; for (i = 0; i < exo->num_nodes; i++) { node_count += node_flag_1ss[i]; } #ifdef PARALLEL k = ProcWithMaxInt(node_count, &i); MPI_Bcast(&(bc_ptr->BC_matrl_index_1), 1, MPI_INT, k, MPI_COMM_WORLD); MPI_Bcast(&(bc_ptr->BC_matrl_index_2), 1, MPI_INT, k, MPI_COMM_WORLD); MPI_Bcast(&(bc_ptr->BC_matrl_index_3), 1, MPI_INT, k, MPI_COMM_WORLD); MPI_Bcast(&(bc_ptr->BC_matrl_index_4), 1, MPI_INT, k, MPI_COMM_WORLD); #ifdef DEBUG_HKM print_sync_start(FALSE); if ( !ProcID ) { printf("Final matrl_index's for ibc = %d:\n", ibc); printf("\tBC_matrl_index_1 = %d\n", bc_ptr->BC_matrl_index_1); printf("\tBC_matrl_index_2 = %d\n", bc_ptr->BC_matrl_index_2); printf("\tBC_matrl_index_3 = %d\n", bc_ptr->BC_matrl_index_3); printf("\tBC_matrl_index_4 = %d\n", bc_ptr->BC_matrl_index_4); fflush(stdout); } print_sync_end(FALSE); #endif #endif } safer_free((void **) &bin_matrl); safer_free((void **) &node_flag_1ss); }
int setup_problem(Exo_DB *exo, /* ptr to the finite element mesh database */ Dpi *dpi) /* distributed processing information */ /******************************************************************** * * setup_problem(): * * Setup_problem() determines the degrees of freedom at each * node and formulates the solution vector. It then determines the * communications pattern for exchanging that solution vector between * processors. * Lastly, it sets up structures that help to carry out the * boundary condition integrals on side sets. * * NOTE: * This function was formed by taking common parts out of **********************************************************************/ { int i; static char *yo = "setup_problem:"; /* * Initialize nodal based structures pertaining to properties * and boundary conditions */ init_nodes(exo, dpi); /* * Determine the side set or node set index that matches all of * those boundary conditions. Store it in the boundary condition * structure. */ bc_set_index(exo); /* * Determine what boundary conditions are on side sets that * are internal-boundary side sets */ bc_internal_boundary(exo); /* * Determine the neighboring material indecises on all boundary * conditions whether or not they be specified in the input deck */ bc_matrl_index(exo); /* * Determine from the boundary conditions and problem description * structures what degrees of freedom are shared between materials * and which degrees of freedom generate multiple degrees of freedom * at each node */ coordinate_discontinuous_variables(exo, dpi); /* * Reconcile the boundary conditions with material properties * entered into the databases. */ reconcile_bc_to_matrl(); /* * Determine the values of the DV_Indexing_Type field for boundary * conditions */ determine_dvi_index(); /* * Enumerate my own degrees of freedom on this processor. */ setup_local_nodal_vars(exo, dpi); #ifdef DEBUG_HKM printf("%s: Proc %d after setup_local_nodal_vars\n", yo, ProcID); fflush(stdout); #endif /* * setup communications patterns between ghost and owned * nodes. */ setup_nodal_comm_map(exo, dpi, cx); #ifdef DEBUG_HKM printf("%s: Proc %d after setup_nodal_comm_map\n", yo, ProcID); fflush(stdout); #endif /* * Find the global maximum number of unknowns located at any one * node on any processor */ MaxVarPerNode = find_MaxUnknownNode(); #ifdef DEBUG_HKM printf("%s: Proc %d after find_MaxUnknownNode\n", yo, ProcID); fflush(stdout); #endif /* * Exchange my idea of what materials I have at each node with * my surrounding processors. Make sure we are all in sync */ setup_external_nodal_matrls(exo, dpi, cx); /* * Exchange my idea of what degrees of freedom I have with my * surrounding processors. Make sure we are all in sync. * Owned nodes tell ghost nodes what variables are active * at that node. */ setup_external_nodal_vars(exo, dpi, cx); #ifdef DEBUG_HKM printf("%s: Proc %d after setup_external_nodal_vars\n", yo, ProcID); fflush(stdout); #endif /* * Finish setting the unknown map on this processor */ set_unknown_map(exo, dpi); /* * Now determine the communications pattern that is necessary to * exchange the solution vector in as efficient a manner as * possible */ #ifdef DEBUG_HKM printf("P_%d: setup_dof_comm_map() begins...\n", ProcID); fflush(stdout); #endif log_msg("setup_dof_comm_map..."); setup_dof_comm_map(exo, dpi, cx); /* * Output some statistics concerning the communications pattern */ if (Num_Proc > 1) output_comm_stats(dpi, cx); #ifdef COUPLED_FILL /* * I extracted this from setup_fill_comm_map because some of the * renormalization routines make use of them. */ num_fill_unknowns = count_vardofs(FILL, dpi->num_universe_nodes); internal_fill_unknowns = count_vardofs(FILL, dpi->num_internal_nodes); owned_fill_unknowns = count_vardofs(FILL, (dpi->num_internal_nodes + dpi->num_boundary_nodes)); boundary_fill_unknowns = owned_fill_unknowns - internal_fill_unknowns; external_fill_unknowns = num_fill_unknowns - owned_fill_unknowns; #else /* COUPLED_FILL */ /* * Set up the communications pattern for doing a solution * of the FILL variable type equations alone, if needed */ if (Explicit_Fill) { #ifdef DEBUG fprintf(stderr, "P_%d: setup_fill_comm_map() begins...\n", ProcID); #endif /* DEBUG */ log_msg("setup_fill_comm_map..."); setup_fill_comm_map(exo, dpi, cx); } #endif /* COUPLED_FILL */ /* * Possibly increase the number of variable descriptions to include * those that are not part of the solution vector */ vdesc_augment(); /* * Setup Boundary condition inter-connectivity * -> Stefan flow bc's need to know what bc's contain the reaction * info. * -> YFLUX bc's need to be connected -> future implementation */ set_up_BC_connectivity(); /* * Set up the structures necessary to carry out * surface integrals */ #ifdef DEBUG fprintf(stderr, "P_%d set_up_Surf_BC() begins...\n",ProcID); #endif log_msg("set_up_Surf_BC..."); set_up_Surf_BC(First_Elem_Side_BC_Array, exo, dpi); /* * Set up the Edge boundary condition structures */ #ifdef DEBUG fprintf(stderr, "P_%d: set_up_Edge_BC() begins...\n",ProcID); #endif log_msg("set_up_Edge_BC..."); set_up_Edge_BC(First_Elem_Edge_BC_Array, exo, dpi); /* * Set up "boundary" conditions on level set surfaces */ set_up_Embedded_BC(); /* * Print out the edge boudary condition structures * if necessary */ if (Debug_Flag) { print_setup_Surf_BC(First_Elem_Side_BC_Array); } /* * Malloc structures of size Num_Var_Info_Records */ ei->VDindex_to_Lvdesc = alloc_int_1(Num_Var_Info_Records, -1); for (i = 0; i < MAX_ELEMENT_INDICES_RELATED; i++) { eiRelated[i]->VDindex_to_Lvdesc = alloc_int_1(Num_Var_Info_Records, -1); } /* * Setup the storage for temporary quantities of interest that * are storred on a "per processor element" basis. These include * temporary storage of volumetric quadrature information */ setup_element_storage(); /* * Setup some structures for solving problems with shell elements. */ init_shell_element_blocks(exo); /* Communicate non-shared but needed BC information */ exchange_bc_info(); return 0; }
int coordinate_discontinuous_variables(Exo_DB *exo, Dpi *dpi) /******************************************************************** * * coordinate_discontinuous_variables(): * * -> Make sure we have the correct designations for the v field * in the problem description structure for each material. * * -> Make sure that we have the same v field for all material * types on all processors. * * *******************************************************************/ { int ibc, eqn_type, ss_index, side_index, k, node_num, imat; int num_mat, mat_index, var_type, *ivec; UMI_LIST_STRUCT *curr_mat_list; NODE_INFO_STRUCT *node_ptr; PROBLEM_DESCRIPTION_STRUCT *curr_pd; /* * Loop over the boundary conditions. If we have a cross * phase discontinuous boundary condition, then we need to set the * v field for the appropriate variable types on both sides of the * interface to denote a discontinuous interpolation at the * interface. */ for (ibc = 0; ibc < Num_BC; ibc++) { if (BC_Types[ibc].desc->i_apply == CROSS_PHASE_DISCONTINUOUS) { eqn_type = BC_Types[ibc].desc->equation; /* * If we are applying a bc on the momentum equations * let's assign it a base equation type */ if (eqn_type == R_MOMENTUM1 || eqn_type == R_MOMENTUM2 || eqn_type == R_MOMENTUM3 || eqn_type == R_MOM_NORMAL || eqn_type == R_MOM_TANG1 || eqn_type == R_MOM_TANG2 ) { eqn_type = R_MOMENTUM1; } /* * If we are applying a discontinuous bc on one species * equation, then we must apply it to all species equations. */ if (eqn_type == R_MASS || (eqn_type >= R_SPECIES_UNK_0 && eqn_type <= R_SPECIES_UNK_LAST) ) { eqn_type = R_SPECIES_UNK_0; } for (ss_index = 0; ss_index < exo->num_side_sets; ss_index++) { if (BC_Types[ibc].BC_ID == exo->ss_id[ss_index]) { for (side_index = 0; side_index < exo->ss_num_sides[ss_index]; side_index++) { for (k = exo->ss_node_side_index[ss_index][side_index]; k < exo->ss_node_side_index[ss_index][side_index+1]; k++) { node_num = exo->ss_node_list[ss_index][k]; node_ptr = Nodes[node_num]; curr_mat_list = &(node_ptr->Mat_List); num_mat = curr_mat_list->Length; /* * Now make sure that we have the discontinuous var turned * on */ for (imat = 0; imat < num_mat; imat++) { mat_index= (curr_mat_list->List)[imat]; curr_pd = pd_glob[mat_index]; if (eqn_type == R_MOMENTUM1) { turn_on_discontinuous(curr_pd, R_MOMENTUM1); turn_on_discontinuous(curr_pd, R_MOMENTUM2); turn_on_discontinuous(curr_pd, R_MOMENTUM3); turn_on_discontinuous(curr_pd, PRESSURE); } else if (eqn_type == R_SPECIES_UNK_0) { turn_on_discontinuous(curr_pd, R_MASS); for (var_type = R_SPECIES_UNK_0; var_type < R_SPECIES_UNK_LAST; var_type++) { turn_on_discontinuous(curr_pd, var_type); } } else { turn_on_discontinuous(curr_pd, eqn_type); } } } } } } } } /* * Just to dot the eyes, make sure that v fields are uniform on * distributed processor problems. We will use the MPI_BOR * operation on a Reduce operation to processor zero, followed * by a broadcast from zero, to accomplish this. */ #ifdef PARALLEL ivec = alloc_int_1(V_LAST, 0); for (imat = 0; imat < upd->Num_Mat; imat++) { curr_pd = pd_glob[imat]; for (k = 0; k < V_LAST; k++) { ivec[k] = curr_pd->v[k]; } ReduceBcast_BOR(ivec, V_LAST); #ifdef DEBUG_HKM print_sync_start(TRUE); for (k = 0; k < V_LAST; k++) { if (curr_pd->v[k] != ivec[k]) { printf("P_%d: v field for var_type %d changed from %d to %d\n", ProcID, k, curr_pd->v[k], ivec[k]); } } print_sync_end(TRUE); #endif for (k = 0; k < V_LAST; k++) { curr_pd->v[k] = ivec[k]; } } safer_free((void **) &ivec); #endif return 0; }
void setup_external_nodal_matrls(Exo_DB *exo, Dpi *dpi, Comm_Ex *cx) /************************************************************************ * * setup_external_nodal_matrls(): * * This routine exchanges information about the materials * from each owned node to each ghost node. We first pack the information * about the solution vector on each owned node on this current * processor into a compact form. Then, we use the normal exchange * node information routine to exchange this information with the * neighboring processors. * Then, we unpack the information obtained from neighboring * processors into nodal variable structures. And, then we use * the same procedure that we used in setup_local_nodal_vars() to * assign nodal var structures to external nodes. * At the end of this procedure, we are assured that ghost nodes will * have the same picture of the solution vector as owned nodes. * ************************************************************************/ { int i, k, p, node_num, max_Matrl; int *mesg_send = NULL, *mesg_recv = NULL, *istart; NODE_INFO_STRUCT *node; COMM_NP_STRUCT *np_ptr, *np_base = NULL; /* * Find out the maximum number of materials per node */ max_Matrl = find_MaxMatrlPerNode(); /* * Pack information for sending */ if (dpi->num_neighbors > 0) { mesg_send = alloc_int_1(ptr_node_send[dpi->num_neighbors] * max_Matrl, -1); mesg_recv = alloc_int_1(ptr_node_recv[dpi->num_neighbors] * max_Matrl, -1); for (i = 0; i < ptr_node_send[dpi->num_neighbors]; i++) { node_num = list_node_send[i]; node = Nodes[node_num]; istart = mesg_send + i*max_Matrl; for (k = 0; k < node->Mat_List.Length; k++) { istart[k] = node->Mat_List.List[k]; } } np_base = alloc_struct_1(COMM_NP_STRUCT, dpi->num_neighbors); } np_ptr = np_base; for (p = 0; p < dpi->num_neighbors; p++) { np_ptr->neighbor_ProcID = cx[p].neighbor_name; np_ptr->send_message_buf = (void *) (mesg_send + ptr_node_send[p] * max_Matrl); np_ptr->send_message_length = sizeof(int) * cx[p].num_nodes_send * max_Matrl; np_ptr->recv_message_buf = (void *) (mesg_recv + ptr_node_recv[p] * max_Matrl); np_ptr->recv_message_length = sizeof(int) * cx[p].num_nodes_recv * max_Matrl; np_ptr++; } exchange_neighbor_proc_info(dpi->num_neighbors, np_base); #ifdef DEBUG_HKM printf("P_%d at barrier after exchange in setup_external_nodal_matrls\n", ProcID); fflush(stdout); #ifdef PARALLEL MPI_Barrier(MPI_COMM_WORLD); #endif #endif for (i = 0; i < dpi->num_external_nodes; i++) { node_num = dpi->num_internal_nodes + dpi->num_boundary_nodes + i; node = Nodes[node_num]; istart = mesg_recv + i*max_Matrl; /* * Unpack ther materials and then add them to the existing * list. */ for (k = 0; k < max_Matrl; k++) { if (istart[k] != -1) { add_to_umi_int_list(&(node->Mat_List), istart[k]); } else { break; } } } #ifdef DEBUG_HKM if ((k = find_MaxMatrlPerNode()) != max_Matrl) { printf("ERROR! %d %d\n", max_Matrl, k); exit (-1); } #endif /* * Free memory allocated in this routine */ safer_free((void **) &np_base); safer_free((void **) &mesg_send); safer_free((void **) &mesg_recv); /* * When in debug mode, print out a complete listing of variables at * every node */ #ifdef DEBUG_HKM printf("P_%d at barrier at end of setup_external_nodal_matrlx\n", ProcID); fflush(stdout); #ifdef PARALLEL MPI_Barrier(MPI_COMM_WORLD); #endif #endif }