/* This routine will perform all of the eigenvalue handling for * "regular" systems. That means it is NOT for 3D stability of a 2D * flow. */ int solve_full_stability_problem(struct Aztec_Linear_Solver_System *ams, double x[], /* Value of the solution vector */ double delta_t, /* Time step size */ double theta, /* Variable time integration parameter * explicit (theta = 1) to implicit * (theta = 0) */ double resid_vector[], double x_old[], /* Value of the old solution vector */ double x_older[], /* */ double xdot[], /* Value of xdot predicted for new * solution */ double xdot_old[], /* */ double x_update[], /* */ int *converged, /* Whether the Newton has converged */ int *nprint, /* Counter for time step number */ int tnv, /* Number of nodal results */ int tnv_post, /* Number of post processing results */ struct Results_Description *rd, int *gindex, int *p_gsize, double *gvec, double time_value, Exo_DB *exo, /* Ptr to finite element mesh db */ Dpi *dpi) /* Ptr to distributed processing info */ { int istuff[PARAMETER_ARRAY_LENGTH]; double dstuff[PARAMETER_ARRAY_LENGTH]; int i, j, mn; int num_total_nodes; /* Number of nodes that each processor is responsible for */ double *scale, *zero; double double_temp; double *jacobian_matrix; double *mass_matrix, *mass_matrix_tmpA, *mass_matrix_tmpB; int *ija = ams->bindx; /* This structure is the same for ALL matrices... */ /* Initialize... */ zero = calloc(NumUnknowns, sizeof(double)); init_vec_value(&zero[0], 0.0, NumUnknowns); scale = calloc(NumUnknowns, sizeof(double)); init_vec_value(scale, 0.0, NumUnknowns); mass_matrix = mass_matrix_tmpA = mass_matrix_tmpB = NULL; jacobian_matrix = NULL; /* Used to let the output routines know that we're not doing normal * mode analysis and just create one set of output files... */ LSA_3D_of_2D_wave_number = -1.0; LSA_current_wave_number = 0; /* Allocate mass matrix */ mass_matrix = calloc((NZeros+5), sizeof(double)); if(LSA_COMPARE) { mass_matrix_tmpA = calloc((NZeros+5), sizeof(double)); mass_matrix_tmpB = calloc((NZeros+5), sizeof(double)); } num_total_nodes = dpi->num_universe_nodes; /* Get jacobian matrix */ printf("Assembling J...\n"); jacobian_matrix = ams->val; init_vec_value(&jacobian_matrix[0], 0.0, (NZeros+1)); af->Assemble_Residual = TRUE; af->Assemble_Jacobian = TRUE; af->Assemble_LSA_Jacobian_Matrix = TRUE; af->Assemble_LSA_Mass_Matrix = FALSE; matrix_fill_full(ams, &x[0], &resid_vector[0], x_old, x_older, xdot, xdot_old, x_update, &delta_t, &theta, First_Elem_Side_BC_Array, &time_value, exo, dpi, &num_total_nodes, &zero[0], &zero[0], NULL); /* This call will fill in scale[] with the proper row scales. I * need to get this in terms of the Jacobian, and the apply it to * the "mass" matrix, later... */ row_sum_scaling_scale(ams, resid_vector, scale); /* Now compute mass matrix. */ /* MMH: * theta = 0 makes the method implicit, so that we get the * var_{}^{n+1} coefficient. */ theta = 0.0; /* The other way: */ if(LSA_COMPARE) /* For a comparison, this is the Subtraction method... */ { init_vec_value(&mass_matrix_tmpA[0], 0.0, (NZeros+1)); init_vec_value(&mass_matrix_tmpB[0], 0.0, (NZeros+1)); TimeIntegration = TRANSIENT; for(mn = 0; mn < upd->Num_Mat; mn++) { pd_glob[mn]->TimeIntegration = TRANSIENT; for(i = 0; i < MAX_EQNS; i++) if(pd_glob[mn]->e[i]) { pd_glob[mn]->e[i] |= T_MASS; pd_glob[mn]->etm[i][LOG2_MASS] = 1.0; } } printf("Assembling B_tmpA...\n"); af->Assemble_Residual = TRUE; af->Assemble_Jacobian = TRUE; af->Assemble_LSA_Jacobian_Matrix = FALSE; af->Assemble_LSA_Mass_Matrix = FALSE; ams->val = mass_matrix_tmpA; delta_t = 1.0; /* To get the phi_i*phi_j terms * (transient soln.) */ tran->delta_t = 1.0; /*for Newmark-Beta terms in Lagrangian Solid*/ matrix_fill_full(ams, x, resid_vector, x_old, x_older, xdot, xdot_old, x_update, &delta_t, &theta, First_Elem_Side_BC_Array, &time_value, exo, dpi, &num_total_nodes, &zero[0], &zero[0], NULL); printf("Assembling B_tmpB...\n"); af->Assemble_Residual = TRUE; af->Assemble_Jacobian = TRUE; af->Assemble_LSA_Jacobian_Matrix = FALSE; af->Assemble_LSA_Mass_Matrix = FALSE; ams->val = mass_matrix_tmpB; delta_t = 1.0e+72; /* To simulate steady-state soln. */ tran->delta_t = 1.0e+72; /*for Newmark-Beta terms in Lagrangian Solid*/ matrix_fill_full(ams, x, resid_vector, x_old, x_older, xdot, xdot_old, x_update, &delta_t, &theta, First_Elem_Side_BC_Array, &time_value, exo, dpi, &num_total_nodes, &zero[0], &zero[0], NULL); for(i = 0; i < NZeros+1; i++) { double_temp = mass_matrix_tmpB[i]-mass_matrix_tmpA[i]; /* if(double_temp != 0.0) if(fabs(double_temp/mass_matrix_tmpA[i]) < 1e-2) printf("Index %d is near-annihilatory. value = %g, difference = %g, ratio = %g\n", i, mass_matrix_tmpA[i], double_temp, fabs(double_temp/mass_matrix_tmpA[i])); */ mass_matrix_tmpA[i] = double_temp; } } /* Initialize */ for(mn = 0; mn < upd->Num_Mat; mn++) { pd_glob[mn]->TimeIntegration = TRANSIENT; for(i = 0; i < MAX_EQNS; i++) /* MMH: Note that there is apparently no support for * multiple species equations. */ if(pd_glob[mn]->e[i]) { /* printf("Equation %d is on.\n",i); */ pd_glob[mn]->e[i] = T_MASS; for (j=0;j<MAX_TERM_TYPES;j++) pd_glob[mn]->etm[i][j] = 0.0; pd_glob[mn]->etm[i][LOG2_MASS] = 1.0; } } /* Hopefully, with all of these settings, this should compute the * mass matrix in mass_matrix and ija */ init_vec_value(&mass_matrix[0], 0.0, (NZeros+1)); printf("Assembling B...\n"); af->Assemble_Residual = TRUE; af->Assemble_Jacobian = TRUE; af->Assemble_LSA_Jacobian_Matrix = FALSE; af->Assemble_LSA_Mass_Matrix = TRUE; ams->val = mass_matrix; delta_t = 1.0; tran->delta_t = 1.0; /*for Newmark-Beta terms in Lagrangian Solid*/ matrix_fill_full(ams, x, resid_vector, x_old, x_older, xdot, xdot_old, x_update, &delta_t, &theta, First_Elem_Side_BC_Array, &time_value, exo, dpi, &num_total_nodes, &zero[0], &zero[0], NULL); /* MMH: Curious... Why NZeros + 1 ??? */ /* for(i = 0; i < NZeros+1; i++) { mass_matrix[i] = -mass_matrix[i]; mass_matrix[i] = mass_matrix_tmpA[i]; } */ /* Now scale the mass matrix with the jacobian scaling we already * computed. */ matrix_scaling(ams, mass_matrix, -1.0, scale); if (LSA_COMPARE) { matrix_scaling(ams, mass_matrix_tmpA, -1.0, scale); compare_mass_matrices(jacobian_matrix, mass_matrix, mass_matrix_tmpA, ija, NumUnknowns); } /* Print jacobian and mass matrices */ if(Linear_Stability == LSA_SAVE || eigen->Eigen_Matrix_Output) output_stability_matrices(mass_matrix, jacobian_matrix, ija, num_total_nodes, NumUnknowns, NZeros); if(Linear_Stability != LSA_SAVE) { /* Call eggroll initiator */ for(i = 0; i < PARAMETER_ARRAY_LENGTH; i++) { istuff[i] = -1; dstuff[i] = -1.0; } eggroll_init(NumUnknowns, NZeros, &istuff[0], &dstuff[0]); /* Call eggroll solver */ eggrollwrap(istuff, dstuff, ija, jacobian_matrix, mass_matrix, x, ExoFileOut, ProblemType, delta_t, theta, x_old, xdot, xdot_old, resid_vector, converged, nprint, tnv, tnv_post, rd, gindex, p_gsize, gvec, time_value, exo, Num_Proc, dpi); } /* Free space, and return ams->val to what it was when we came * in... */ free(mass_matrix); if(LSA_COMPARE) { free(mass_matrix_tmpA); free(mass_matrix_tmpB); } free(scale); free(zero); ams->val = jacobian_matrix; return(0); }
/* This routine will perform all of the eigenvalue handling for * calculations of 3D stability of a 2D flow. */ int solve_3D_of_2D_stability_problem(struct Aztec_Linear_Solver_System *ams, double x[], /* Value of the solution vector */ double delta_t, /* Time step size */ double theta, /* Variable time integration parameter * explicit (theta = 1) to implicit * (theta = 0) */ double resid_vector[], double x_old[], /* Value of the old solution vector */ double x_older[], /* */ double xdot[], /* Value of xdot predicted for new * solution */ double xdot_old[], /* */ double x_update[], /* */ int *converged, /* Whether the Newton has converged */ int *nprint, /* Counter for time step number */ int tnv, /* Number of nodal results */ int tnv_post, /* Number of post processing results */ struct Results_Description *rd, int *gindex, int *p_gsize, double *gvec, double time_value, Exo_DB *exo, /* Ptr to finite element mesh db */ Dpi *dpi) /* Ptr to distributed processing info */ { int istuff[PARAMETER_ARRAY_LENGTH]; double dstuff[PARAMETER_ARRAY_LENGTH]; int i, j, mn, wn; int num_total_nodes; /* Number of nodes that each processor is responsible for */ double *scale, *zero; double *jacobian_matrix, *jacobian_matrix_1, *jacobian_matrix_2; double *mass_matrix, *mass_matrix_1, *mass_matrix_2, *mass_matrix_1_tmpA, *mass_matrix_1_tmpB, *mass_matrix_2_tmpA, *mass_matrix_2_tmpB; int *ija = ams->bindx; /* This structure is the same for ALL matrices... */ int **e_save; dbl ***etm_save; /* Initialize... */ zero = calloc(NumUnknowns, sizeof(double)); init_vec_value(&zero[0], 0.0, NumUnknowns); scale = calloc(NumUnknowns, sizeof(double)); init_vec_value(scale, 0.0, NumUnknowns); mass_matrix = mass_matrix_1 = mass_matrix_2 = mass_matrix_1_tmpA = mass_matrix_1_tmpB = mass_matrix_2_tmpA = mass_matrix_2_tmpB = NULL; jacobian_matrix = jacobian_matrix_1 = jacobian_matrix_2 = NULL; e_save = (int **)array_alloc(2, upd->Num_Mat, MAX_EQNS, sizeof(int)); etm_save = (dbl ***)array_alloc(3, upd->Num_Mat, MAX_EQNS, MAX_TERM_TYPES, sizeof(dbl)); for(mn = 0; mn < upd->Num_Mat; mn++) for(i = 0; i < MAX_EQNS; i++) { e_save[mn][i] = pd_glob[mn]->e[i]; for(j = 0; j < MAX_TERM_TYPES; j++) etm_save[mn][i][j] = pd_glob[mn]->etm[i][j]; } /* Allocate mass matrices */ mass_matrix = calloc((NZeros+5), sizeof(double)); mass_matrix_1 = calloc((NZeros+5), sizeof(double)); mass_matrix_2 = calloc((NZeros+5), sizeof(double)); if(LSA_COMPARE) { mass_matrix_1_tmpA = calloc((NZeros+5), sizeof(double)); mass_matrix_1_tmpB = calloc((NZeros+5), sizeof(double)); mass_matrix_2_tmpA = calloc((NZeros+5), sizeof(double)); mass_matrix_2_tmpB = calloc((NZeros+5), sizeof(double)); } /* Allocate jacobian matrices */ jacobian_matrix = calloc((NZeros+5), sizeof(double)); /* jacobian_matrix_1 is going to use the ams->val space already * allocated... */ jacobian_matrix_1 = ams->val; jacobian_matrix_2 = calloc((NZeros+5), sizeof(double)); num_total_nodes = dpi->num_universe_nodes; for(wn = 0; wn < LSA_number_wave_numbers; wn++) { LSA_current_wave_number = wn; LSA_3D_of_2D_wave_number = LSA_wave_numbers[wn]; printf("Solving for wave number %d out of %d. WAVE NUMBER = %g\n", wn + 1, LSA_number_wave_numbers, LSA_3D_of_2D_wave_number); /* Get the original pd_glob[mn]->e[i] and pd_glob[mn]->etm[i][*] * values back. */ TimeIntegration = STEADY; theta = 0.0; for(mn = 0; mn < upd->Num_Mat; mn++) for(i = 0; i < MAX_EQNS; i++) { pd_glob[mn]->TimeIntegration = STEADY; pd_glob[mn]->e[i] = e_save[mn][i]; for(j = 0; j < MAX_TERM_TYPES; j++) pd_glob[mn]->etm[i][j] = etm_save[mn][i][j]; } /* Jacobian matrix, pass 1 */ printf("Assembling J (pass 1)...\n"); LSA_3D_of_2D_pass = 1; ams->val = jacobian_matrix_1; init_vec_value(&jacobian_matrix_1[0], 0.0, (NZeros+1)); af->Assemble_Residual = TRUE; af->Assemble_Jacobian = TRUE; af->Assemble_LSA_Jacobian_Matrix = TRUE; af->Assemble_LSA_Mass_Matrix = FALSE; matrix_fill_full(ams, x, resid_vector, x_old, x_older, xdot, xdot_old, x_update, &delta_t, &theta, First_Elem_Side_BC_Array, &time_value, exo, dpi, &num_total_nodes, &zero[0], &zero[0], NULL); /* Jacobian matrix, pass 2 */ printf("Assembling J (pass 2)...\n"); LSA_3D_of_2D_pass = 2; ams->val = jacobian_matrix_2; init_vec_value(&jacobian_matrix_2[0], 0.0, (NZeros+1)); af->Assemble_Residual = TRUE; af->Assemble_Jacobian = TRUE; af->Assemble_LSA_Jacobian_Matrix = TRUE; af->Assemble_LSA_Mass_Matrix = FALSE; matrix_fill_full(ams, x, resid_vector, x_old, x_older, xdot, xdot_old, x_update, &delta_t, &theta, First_Elem_Side_BC_Array, &time_value, exo, dpi, &num_total_nodes, &zero[0], &zero[0], NULL); /* Add the two passes together to get the "real" Jacobian * matrix. */ for(i = 0; i < NZeros+1; i++) jacobian_matrix[i] = jacobian_matrix_1[i] + jacobian_matrix_2[i]; /* This call will fill in scale[] with the proper row scales. I * need to get this in terms of the Jacobian, and the apply it to * the "mass" matrix, later... */ ams->val = jacobian_matrix; row_sum_scaling_scale(ams, resid_vector, scale); /* Now compute the mass matrix. This is more complicated than * it seems if we want to perform a comparison with the * Subtractionism method. There are, in fact, three different * matrices with which we could compare along the way: mass * matrix from pass 1, mass matrix from pass 2, and final * (Additionism?) mass matrix. I will construct comparisons * for all three... * * There is an order restriction. We cannot recover the * pd_glob[mn]->e[i] (and pd_glob[mn]->etm[i][*]) values AFTER * we compute the mass matrix with my method, because we set * everything to 0 except the T_MASS-ish terms. So we should * compute all of the Subtractionism matrices, * mass_matrix_[1,2]_tmp[A,b], first. Then compute the two * mass_matrix_[1,2] passes. Since we save the pd_glob[mn] * structures to reset things over multiple wave number * computations, we could refer to the saved values to remove * this restriction, but there's less code this way... */ /* theta = 0 makes the method implicit, so that we get the * var_{}^{n+1} coefficient. We want this for ALL mass matrix * computations, but not for the jacobian computations (they * should have the mass etm = 0, so it wouldn't matter...). */ theta = 0.0; TimeIntegration = TRANSIENT; for(mn = 0; mn < upd->Num_Mat; mn++) { pd_glob[mn]->TimeIntegration = TRANSIENT; for(i = 0; i < MAX_EQNS; i++) if(pd_glob[mn]->e[i]) { pd_glob[mn]->e[i] |= T_MASS; pd_glob[mn]->etm[i][LOG2_MASS] = 1.0; } } /* Subtractionism method for mass matrix, passes 1 and 2. */ if(LSA_COMPARE) { init_vec_value(&mass_matrix_1_tmpA[0], 0.0, (NZeros+1)); init_vec_value(&mass_matrix_1_tmpB[0], 0.0, (NZeros+1)); init_vec_value(&mass_matrix_2_tmpA[0], 0.0, (NZeros+1)); init_vec_value(&mass_matrix_2_tmpB[0], 0.0, (NZeros+1)); /* They all have the same action flag settings... */ af->Assemble_Residual = TRUE; af->Assemble_Jacobian = TRUE; af->Assemble_LSA_Jacobian_Matrix = FALSE; af->Assemble_LSA_Mass_Matrix = FALSE; printf("Assembling B_tmpA (pass 1)...\n"); LSA_3D_of_2D_pass = 1; ams->val = mass_matrix_1_tmpA; delta_t = 1.0; /* To get the phi_i*phi_j terms * (transient soln.) */ tran->delta_t = 1.0; /*for Newmark-Beta terms in Lagrangian Solid*/ matrix_fill_full(ams, x, resid_vector, x_old, x_older, xdot, xdot_old, x_update, &delta_t, &theta, First_Elem_Side_BC_Array, &time_value, exo, dpi, &num_total_nodes, &zero[0], &zero[0], NULL); printf("Assembling B_tmpB (pass 1)...\n"); LSA_3D_of_2D_pass = 1; ams->val = mass_matrix_1_tmpB; delta_t = 1.0e+72; /* To simulate steady-state soln. */ tran->delta_t = 1.0e+72; /*for Newmark-Beta terms in Lagrangian Solid*/ matrix_fill_full(ams, x, resid_vector, x_old, x_older, xdot, xdot_old, x_update, &delta_t, &theta, First_Elem_Side_BC_Array, &time_value, exo, dpi, &num_total_nodes, &zero[0], &zero[0], NULL); printf("Assembling B_tmpA (pass 2)...\n"); LSA_3D_of_2D_pass = 2; ams->val = mass_matrix_2_tmpA; delta_t = 1.0; /* To get the phi_i*phi_j terms * (transient soln.) */ tran->delta_t = 1.0; /*for Newmark-Beta terms in Lagrangian Solid*/ matrix_fill_full(ams, x, resid_vector, x_old, x_older, xdot, xdot_old, x_update, &delta_t, &theta, First_Elem_Side_BC_Array, &time_value, exo, dpi, &num_total_nodes, &zero[0], &zero[0], NULL); printf("Assembling B_tmpB (pass 2)...\n"); LSA_3D_of_2D_pass = 2; ams->val = mass_matrix_2_tmpB; delta_t = 1.0e+72; /* To simulate steady-state soln. */ tran->delta_t = 1.0e+72; /*for Newmark-Beta terms in Lagrangian Solid*/ matrix_fill_full(ams, &x[0], &resid_vector[0], x_old, x_older, xdot, xdot_old, x_update, &delta_t, &theta, First_Elem_Side_BC_Array, &time_value, exo, dpi, &num_total_nodes, &zero[0], &zero[0], NULL); /* Now construct the actual Subtractionism versions of each * pass. These are stored in mass_matrix_[1,2]_tmpA. */ for(i = 0; i < NZeros+1; i++) { mass_matrix_1_tmpA[i] = mass_matrix_1_tmpB[i] - mass_matrix_1_tmpA[i]; mass_matrix_2_tmpA[i] = mass_matrix_2_tmpB[i] - mass_matrix_2_tmpA[i]; } } /* Now we want to compute both passes of our mass matrix using * the action flag. We already set some of the parameters * correcetly above... */ /* Initialize */ for(mn = 0; mn < upd->Num_Mat; mn++) { for(i = 0; i < MAX_EQNS; i++) if(pd_glob[mn]->e[i]) { pd_glob[mn]->e[i] = T_MASS; for (j = 0; j < MAX_TERM_TYPES; j++) pd_glob[mn]->etm[i][j] = 0.0; pd_glob[mn]->etm[i][LOG2_MASS] = 1.0; } } init_vec_value(&mass_matrix[0], 0.0, (NZeros+1)); init_vec_value(&mass_matrix_1[0], 0.0, (NZeros+1)); init_vec_value(&mass_matrix_2[0], 0.0, (NZeros+1)); /* These action flags and delta_t are the same for both * passes. */ af->Assemble_Residual = TRUE; af->Assemble_Jacobian = TRUE; af->Assemble_LSA_Jacobian_Matrix = FALSE; af->Assemble_LSA_Mass_Matrix = TRUE; delta_t = 1.0; tran->delta_t = 1.0; /*for Newmark-Beta terms in Lagrangian Solid*/ printf("Assembling B (pass 1)...\n"); LSA_3D_of_2D_pass = 1; ams->val = mass_matrix_1; matrix_fill_full(ams, x, resid_vector, x_old, x_older, xdot, xdot_old, x_update, &delta_t, &theta, First_Elem_Side_BC_Array, &time_value, exo, dpi, &num_total_nodes, &zero[0], &zero[0], NULL); printf("Assembling B (pass 2)...\n"); LSA_3D_of_2D_pass = 2; ams->val = mass_matrix_2; matrix_fill_full(ams, x, resid_vector, x_old, x_older, xdot, xdot_old, x_update, &delta_t, &theta, First_Elem_Side_BC_Array, &time_value, exo, dpi, &num_total_nodes, &zero[0], &zero[0], NULL); /* Get my mass matrix by adding the two passes together (don't * forget the -! */ for(i = 0; i < NZeros+1; i++) mass_matrix[i] = -(mass_matrix_1[i] + mass_matrix_2[i]); /* Scale all of the preliminary pieces for comparison, if we're * doing one. Otherwise, just scale the above mass_matrix... */ if(LSA_COMPARE) { for(i = 0; i < NumUnknowns; i++) { mass_matrix_1_tmpA[i] /= scale[i]; mass_matrix_2_tmpA[i] /= scale[i]; mass_matrix_1[i] /= scale[i]; mass_matrix_2[i] /= scale[i]; mass_matrix[i] /= scale[i]; for(j = ija[i]; j < ija[i+1]; j++) { mass_matrix_1_tmpA[j] /= scale[i]; mass_matrix_2_tmpA[j] /= scale[i]; mass_matrix_1[j] /= scale[i]; mass_matrix_2[j] /= scale[i]; mass_matrix[j] /= scale[i]; } } /* Compare 1st pass results. */ printf("Comparing pass 1 of mass matrix with Subtractionism method..."); compare_mass_matrices(jacobian_matrix, mass_matrix_1, mass_matrix_1_tmpA, ija, NumUnknowns); /* Comapre 2nd pass results. */ printf("Comparing pass 2 of mass matrix with Subtractionism method..."); compare_mass_matrices(jacobian_matrix, mass_matrix_2, mass_matrix_2_tmpA, ija, NumUnknowns); /* Get the sum of passes for Subtractionism method. This * was already done above for my regular method. */ for(i = 0; i < NZeros+1; i++) mass_matrix_1_tmpA[i] = mass_matrix_1_tmpA[i] + mass_matrix_1_tmpB[i]; /* Compare sum of passes results. */ printf("Comparing sum of passes of mass matrix with Subtractionism method..."); compare_mass_matrices(jacobian_matrix, mass_matrix, mass_matrix_1_tmpA, ija, NumUnknowns); } else for(i = 0; i < NumUnknowns; i++) { mass_matrix[i] /= scale[i]; for(j = ija[i]; j < ija[i+1]; j++) mass_matrix[j] /= scale[i]; } /* Print Jacobian and mass matrices */ if(Linear_Stability == LSA_3D_OF_2D_SAVE || eigen->Eigen_Matrix_Output) output_stability_matrices(mass_matrix, jacobian_matrix, ija, num_total_nodes, NumUnknowns, NZeros); if(Linear_Stability != LSA_3D_OF_2D_SAVE) { /* Call eggroll initiator */ for(i = 0; i < PARAMETER_ARRAY_LENGTH; i++) { istuff[i] = -1; dstuff[i] = -1.0; } printf("Entering eggrollinit...\n"); fflush(stdout); eggroll_init(NumUnknowns, NZeros, &istuff[0], &dstuff[0]); /* Call eggroll solver */ printf("Entering eggrollwrap...\n"); fflush(stdout); eggrollwrap(istuff, dstuff, ija, jacobian_matrix, mass_matrix, x, ExoFileOut, ProblemType, delta_t, theta, x_old, xdot, xdot_old, resid_vector, converged, nprint, tnv, tnv_post, rd, gindex, p_gsize, gvec, time_value, exo, Num_Proc, dpi); } } /* for(wn = 0; wn < LSA_num_wave_numbers; wn++) */ /* Free space, and return ams->val to what it was when we came * in... */ free(jacobian_matrix_2); free(jacobian_matrix); if(LSA_COMPARE) { free(mass_matrix_2_tmpB); free(mass_matrix_2_tmpA); free(mass_matrix_1_tmpB); free(mass_matrix_1_tmpA); } free(mass_matrix_2); free(mass_matrix_1); free(mass_matrix); free(e_save); free(etm_save); free(scale); free(zero); ams->val = jacobian_matrix_1; return(0); }
void extract_nodal_vec(double sol_vec[], int var_no, int ktype, int matIndex, double nodal_vec[], Exo_DB *exo, int timeDerivative, double time) /**************************************************************************** * * This function puts the nodal values of the selected variable * into a global solution vector which contains all the mesh nodes, * and interpolates the mid-side and centroid values of * all variables with Q1 interpolation on a * 9-NODE mesh and interpolates to find their values at the mid-side * nodes and at the centroid. * * Now this is set up to be at least compatible w/ parallel computing. * We actually load the nodal vector for the current processor only. * * Written by: Rich Cairncross 27 July 1995 * * Revised: 1997/08/26 14:52 MDT [email protected] * * Input * -------- * sol_vec[] = Current global solution vector * var_no = Variable type to be extracted * ktype = sub_var number of the variable type to be extracted * matIndex = material index of the variable to be extracted. * -1 : extract the nonspecific variable * -2 : extract the first variable with var_no no * matter what material index. * exo = Exodus database structure * timeDerivative = Are we extracting a time derivative? if so, * then this is true. If not, false. * * Output * ------- * nodal_vec[] = nodal vector which receives the value * of the extracted vector. (length number * of nodes) * * * Special Case Processing: * * If var_no is a Species variable, and if we are * solving for nondilute mass transfer, and thus the * number of species may not be equal to the number of species * equations, then this routine can extract the implicit * value of the species variable dropped from the equation * set. * nodal_vec[] in this case for this unknown will be equal * to sum of all of the other species unknowns in the problem. *************************************************************************/ { int eb_index, mn; /* * Initialize the output vector to zero. Thus, we don't have to * visit each node to get a valid result */ init_vec_value(nodal_vec, 0.0, Num_Node); /* * Loop backwards through the element blocks */ for (eb_index = exo->num_elem_blocks - 1; eb_index >= 0; eb_index--) { mn = Matilda[eb_index]; pd = pd_glob[mn]; /* * First check to see whether the variable is defined to be * in the solution vector in the element block * -> If it isn't defined, we don't need to go into the eb routine, * and we will avoid some errors with cross-phase situations. */ if (pd->e[var_no]) { /* * Next check to see whether we are obtaining material * specific values or general values of the variable */ if (matIndex < 0 || matIndex == mn) { /* * Now go into the element block and get the variable at * the nodes. */ extract_nodal_eb_vec(sol_vec, var_no, ktype, matIndex, eb_index, nodal_vec, exo, timeDerivative, time); } } } return; }