/* ********************************************************************* */ unsigned char CheckZone (int z, int bit) /*! * Check if bit "bit" is turned on in zone "z". * *********************************************************************** */ { if (flag == NULL){ /* ---------------------------------------------- the flag array is not defined when Chombo is converting conservative variables to primitive ones. ---------------------------------------------- */ #ifdef CH_SPACEDIM return (0); #endif print1 ("! CheckZone: NULL pointer in CheckZone\n"); QUIT_PLUTO(1); } if (g_dir == IDIR) return (flag[g_k][g_j][z] & bit); if (g_dir == JDIR) return (flag[g_k][z][g_i] & bit); if (g_dir == KDIR) return (flag[z][g_j][g_i] & bit); else{ print1 ("! CheckZone: wrong direction\n"); QUIT_PLUTO(1); } }
/* ********************************************************************* */ void GetOutputFrequency(Output *output, char *output_format) /*! * Set the intervals between output files. * This can be done in three different ways: * * - dt: time interval in code units * - dn: step interval * - dclock: actual clock time (in hours) * * However, dn and dclock are mutually exclusive. * *********************************************************************** */ { char *str; int len, nhrs, nmin, nsec; /* -- time interval in code units (dt) -- */ output->dt = atof(ParamFileGet(output_format, 1)); /* -- check the 2nd field and decide to set "dn" or "dclock" -- */ str = ParamFileGet(output_format,2); len = strlen(str); if (str[len-1] == 'h'){ output->dclock = atof(str); /* clock interval in hours */ nhrs = (int)output->dclock; /* integer part */ nmin = (int)((output->dclock - nhrs)*100.0); /* remainder in minutes */ if (nmin >= 60){ printf ("! OutputFrequency: number of minutes exceeds 60 in %s output\n", output_format); QUIT_PLUTO(1); } output->dclock = nhrs*3600.0 + nmin*60; /* convert to seconds */ output->dn = -1; }else if (str[len-1] == 'm'){ output->dclock = atof(str); /* clock interval in minutes */ nmin = (int)output->dclock; /* integer part */ nsec = (int)((output->dclock - nmin)*100.0); /* remainder in seconds */ if (nsec >= 60){ printf ("! OutputFrequency: number of seconds exceeds 60 in %s output\n", output_format); QUIT_PLUTO(1); } output->dclock = nmin*60.0 + nsec; output->dn = -1; }else if (str[len-1] == 's'){ output->dclock = atof(str); /* clock interval in seconds */ output->dn = -1; }else{ output->dclock = -1.0; output->dn = atoi(ParamFileGet(output_format, 2)); } }
/* ********************************************************************* */ void SoundSpeed2 (double **v, double *cs2, double *h, int beg, int end, int pos, Grid *grid) /*! * Define the square of the sound speed. * * \param [in] v 1D array of primitive quantities * \param [out] cs2 1D array containing the square of the sound speed * \param [in] h 1D array of enthalpy values * \param [in] beg initial index of computation * \param [in] end final index of computation * \param [in] pos an integer specifying the spatial position * inside the cell (only for spatially-dependent EOS) * \param [in] grid pointer to an array of Grid structures * * \return This function has no return value. *********************************************************************** */ { int i; double theta; #if PHYSICS == RHD || PHYSICS == RMHD Enthalpy(v, h, beg, end); for (i = beg; i <= end; i++) { theta = v[i][PRS]/v[i][RHO]; #if EOS == IDEAL cs2[i] = g_gamma*theta/h[i]; #elif EOS == TAUB cs2[i] = theta/(3.0*h[i])*(5.0*h[i] - 8.0*theta)/(h[i] - theta); #endif } #else print ("! SoundSpeed2: Taub EOS not defined for this physics module.\n"); QUIT_PLUTO(1); #endif }
/* ********************************************************************* */ void Entropy (double **v, double *s, int is, int ie) /*! * Compute the entropy. * * \param [in] v 1D array of primitive quantities * \param [in] s 1D array of entropy values * \param [in] is initial index of computation * \param [in] ie final index of computation * * \return This function has no return value. *********************************************************************** */ { int i; double rho; #if EOS == IDEAL for (i = is; i <= ie; i++){ rho = v[i][RHO]; s[i] = v[i][PRS]/pow(rho,g_gamma); } #elif EOS == ISOTHERMAL || EOS == BAROTROPIC print (" Entropy not defined in isothermal or barotropic MHD\n"); QUIT_PLUTO(1); #endif }
/* ********************************************************************* */ void SoundSpeed2 (double **u, double *cs2, double *h, int beg, int end, int pos, Grid *grid) /*! * Define the square of the sound speed. * * \param [in] u 1D array of primitive quantities * \param [out] cs2 1D array containing the square of the sound speed * \param [in] h 1D array of enthalpy values * \param [in] beg initial index of computation * \param [in] end final index of computation * \param [in] pos an integer specifying the spatial position * inside the cell (only for spatially-dependent EOS) * \param [in] grid pointer to an array of Grid structures * * \return This function has no return value. *********************************************************************** */ { int i; #if EOS == IDEAL for (i = beg; i <= end; i++) cs2[i] = g_gamma*u[i][PRS]/u[i][RHO]; #elif EOS == ISOTHERMAL { int j,k; /* -- used as multidimensional indices -- */ double *x1, *x2, *x3; x1 = grid[IDIR].x; x2 = grid[JDIR].x; x3 = grid[KDIR].x; i = *g_i; j = *g_j; k = *g_k; if (g_dir == IDIR) { x1 = (pos == FACE_CENTER ? grid[IDIR].xr : grid[IDIR].x); for (i = beg; i <= end; i++) cs2[i] = g_isoSoundSpeed*g_isoSoundSpeed; }else if (g_dir == JDIR){ x2 = (pos == FACE_CENTER ? grid[JDIR].xr : grid[JDIR].x); for (j = beg; j <= end; j++) cs2[j] = g_isoSoundSpeed*g_isoSoundSpeed; }else if (g_dir == KDIR){ x3 = (pos == FACE_CENTER ? grid[KDIR].xr : grid[KDIR].x); for (k = beg; k <= end; k++) cs2[k] = g_isoSoundSpeed*g_isoSoundSpeed; } } #elif EOS == BAROTROPIC for (i = is ; i <= ie ; i++) { cs2[i] = g_gamma*BAROTROPIC_PR(u[i][RHO])/u[i][RHO]; /* true only if p = K rho^\gamma */ } #else print ("! SoundSpeed2: not defined for this EoS\n"); QUIT_PLUTO(1); #endif }
/* ********************************************************************* */ void Enthalpy (double **v, real *h, int beg, int end) /*! * Compute the enthalpy. * * \param [in] v 1D array of primitive quantities * \param [in] h 1D array of enthalpy values * \param [in] beg initial index of computation * \param [in] end final index of computation * * \return This function has no return value. *********************************************************************** */ { int i; double gmmr; #if EOS == IDEAL gmmr = g_gamma/(g_gamma - 1.0); for (i = beg; i <= end; i++){ h[i] = gmmr*v[i][PRS]/v[i][RHO]; } #elif EOS == ISOTHERMAL || EOS == BAROTROPIC print (" Enthalpy not defined in isothermal or barotropic MHD\n"); QUIT_PLUTO(1); #endif }
/* ***************************************************************** */ void Radiat (double *v, double *rhs) /*! * Provide r.h.s. for tabulated cooling. * ******************************************************************* */ { int klo, khi, kmid; static int ntab; double mu, T, Tmid, scrh, dT, prs; static double *L_tab, *T_tab, E_cost; FILE *fcool; /* ------------------------------------------- Read tabulated cooling function ------------------------------------------- */ if (T_tab == NULL){ print1 (" > Reading table from disk...\n"); fcool = fopen("cooltable.dat","r"); if (fcool == NULL){ print1 ("! Radiat: cooltable.dat could not be found.\n"); QUIT_PLUTO(1); } L_tab = ARRAY_1D(20000, double); T_tab = ARRAY_1D(20000, double); ntab = 0; while (fscanf(fcool, "%lf %lf\n", T_tab + ntab, L_tab + ntab)!=EOF) { ntab++; } E_cost = UNIT_LENGTH/UNIT_DENSITY/pow(UNIT_VELOCITY, 3.0); }
/* ********************************************************************* */ void Enthalpy (double **uprim, double *h, int beg, int end) /*! * Compute the enthalpy. * * \param [in] v 1D array of primitive quantities * \param [in] h 1D array of enthalpy values * \param [in] beg initial index of computation * \param [in] end final index of computation * * \return This function has no return value. *********************************************************************** */ { int i; double gmmr; #if EOS == IDEAL gmmr = g_gamma/(g_gamma - 1.0); for (i = beg; i <= end; i++){ h[i] = gmmr*uprim[i][PRS]/uprim[i][RHO]; } #elif EOS == ISOTHERMAL print ("! Enthalpy: not defined for isothermal EoS\n"); QUIT_PLUTO(1); #endif }
/* *********************************************************************** */ void PrimRHS (double *v, double *dv, double cs2, double h, double *Adv) /*! * Compute the matrix-vector multiplication \f$ A(\mathbf{v})\cdot * d\mathbf{v} \f$ where A is the matrix of the quasi-linear form * of the MHD equations. * * \b References * * - "A solution adaptive upwind scheme for ideal MHD", * Powell et al., JCP (1999) 154, 284 * * - "An unsplit Godunov method for ideal MHD via constrained transport" * Gardiner \& Stone, JCP (2005) 205, 509 * * \param [in] v vector of primitive variables * \param [in] dv limited (linear) slopes * \param [in] cs2 local sound speed * \param [in] h local enthalpy * \param [out] AdV matrix-vector product * * \return This function has no return value. ************************************************************************* */ { int nv; double tau, scrh; double ch2; tau = 1.0/v[RHO]; /* --------------------------------------------- Adv[k] Contains A[k][*]*dv[*] --------------------------------------------- */ Adv[RHO] = v[VXn]*dv[RHO] + v[RHO]*dv[VXn]; scrh = EXPAND(0.0, + v[BXt]*dv[BXt], + v[BXb]*dv[BXb]); #if EOS == IDEAL Adv[VXn] = v[VXn]*dv[VXn] + tau*(dv[PRS] + scrh); #elif EOS == BAROTROPIC Adv[VXn] = v[VXn]*dv[VXn] + tau*(cs2*dv[RHO] + scrh); #elif EOS == ISOTHERMAL Adv[VXn] = v[VXn]*dv[VXn] + tau*(cs2*dv[RHO] + scrh); #else print ("! PRIM_RHS: not defined for this EoS\n"); QUIT_PLUTO(1); #endif EXPAND( ; ,
/* ********************************************************************* */ void ParabolicFlux (Data_Arr V, const State_1D *state, double **dcoeff, int beg, int end, Grid *grid) /*! * Add the diffusion fluxes to the upwind fluxes for explicit time * integration. * * \param [in] V pointer to the 3D array of cell-centered primitive * variables * \param [in,out] state pointer to a State_1D structure * \param [out] dcoeff the diffusion coefficients 1D array * \param [in] beg initial index of computation * \param [in] end final index of computation * \param [in] grid pointer to an array of Grid structures * *********************************************************************** */ { int i, j, k, nv; static double *par_Eflx; /* ------------------------------------------------------------ define par_Eflx as the sum of all parabolic contributions (i.e., visc+tc+resistive) to the total energy flux. ------------------------------------------------------------ */ if (par_Eflx == NULL) par_Eflx = ARRAY_1D(NMAX_POINT, double); for (i = beg; i <= end + 1; i++){ /* -- use memset instead -- */ for (nv = NVAR; nv--; ) { state->par_src[i][nv] = 0.0; state->par_flx[i][nv] = 0.0; } par_Eflx[i] = 0.0; } /* ------------------------------------------------- 1. Viscosity ------------------------------------------------- */ #if VISCOSITY == EXPLICIT #ifdef FARGO print ("! ParabolicFlux: FARGO incompatible with explicit viscosity.\n"); print (" Try STS or RKC instead\n"); QUIT_PLUTO(1); #endif ViscousFlux (V, state->par_flx, state->par_src, dcoeff, beg, end, grid); for (i = beg; i <= end; i++){ EXPAND(state->flux[i][MX1] += state->par_flx[i][MX1]; , state->flux[i][MX2] += state->par_flx[i][MX2]; ,
/* ********************************************************************* */ int CheckNaN (double **u, int is, int ie, int id) /*! * Cheeck whether the array \c u contains Not-a-Number * (NaN) * *********************************************************************** */ { int ii, nv, i, j; for (ii = is; ii <= ie; ii++) { for (nv = 0; nv < NVAR; nv++) { if (u[ii][nv] != u[ii][nv]) { print (" > NaN found (%d), |", id); Show (u, ii); QUIT_PLUTO(1); } }} return (0); }
/* *************************************************************** */ void Enthalpy (real **uprim, real *h, int beg, int end) /* * * * ***************************************************************** */ { int i; double g_gammar; #if EOS == IDEAL g_gammar = g_gamma/(g_gamma - 1.0); for (i = beg; i <= end; i++){ h[i] = g_gammar*uprim[i][PRS]/uprim[i][RHO]; } #elif EOS == ISOTHERMAL print (" Enthalpy not defined for isothermal EoS\n"); QUIT_PLUTO(1); #endif }
/* *************************************************************** */ void ENTROPY (real **v, real *s, int is, int ie) /* * * * ***************************************************************** */ { int i; double rho; #if EOS == IDEAL for (i = is; i <= ie; i++){ rho = v[i][RHO]; s[i] = v[i][PRS]/pow(rho,g_gamma); } #elif EOS == ISOTHERMAL || EOS == BAROTROPIC print (" Entropy not defined in isothermal or barotropic MHD\n"); QUIT_PLUTO(1); #endif }
/* ********************************************************************* */ void HancockStep (const State_1D *state, int beg, int end, Grid *grid) /*! * Use Taylor expansion to compute time-centered states in primitive * variables. * * \param [in,out] state pointer to a State_1D structure * \param [in] beg initial index of computation * \param [in] end final index of computation * \param [in grid pointer to array of Grid structures * * \b References * - "Riemann Solvers and Numerical Methods for Fluid Dynamics" \n * E.F. Toro * *********************************************************************** */ { int nv, i; double scrh, dt_2, ch2; double Adv[NVAR], dv[NVAR]; double *vp, *vm, *vc; static double **src, *d_dl; /* ----------------------------------------- Check scheme compatibility ----------------------------------------- */ #if INTERPOLATION != LINEAR print1 ("! MUSCL-Hancock scheme works with Linear reconstruction only\n"); QUIT_PLUTO(1); #endif if (src == NULL){ src = ARRAY_2D(NMAX_POINT, NVAR, double); d_dl = ARRAY_1D(NMAX_POINT, double); }
/* ********************************************************************* */ double NextTimeStep (Time_Step *Dts, struct INPUT *ini, Grid *grid) /*! * Compute and return the time step for the next time level * using the information from the previous integration * (Dts->inv_dta and Dts->inv_dp). * * \param [in] Dts pointer to the Time_Step structure * \param [in] ini pointer to the Input structure * \param [in] grid pointer to array of Grid structures * * \return The time step for next time level *********************************************************************** */ { int idim; double dt_adv, dt_par, dtnext; double dxmin; double xloc, xglob; /* --------------------------------------------------- 1. Take the maximum of inv_dt across all processors --------------------------------------------------- */ #ifdef PARALLEL xloc = Dts->inv_dta; MPI_Allreduce (&xloc, &xglob, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD); Dts->inv_dta = xglob; #if (PARABOLIC_FLUX != NO) xloc = Dts->inv_dtp; MPI_Allreduce (&xloc, &xglob, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD); Dts->inv_dtp = xglob; #endif #if COOLING != NO xloc = Dts->dt_cool; MPI_Allreduce (&xloc, &xglob, 1, MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD); Dts->dt_cool = xglob; #endif #endif /* ---------------------------------- 2. Compute time step ---------------------------------- */ #if (PARABOLIC_FLUX & EXPLICIT) dt_adv = 1.0/(Dts->inv_dta + 2.0*Dts->inv_dtp); #else dt_adv = 1.0/Dts->inv_dta; #endif dt_adv *= ini->cfl; dtnext = dt_adv; /* ------------------------------------------------------- 3. Maximum propagation speed for the local processor. Global glm_ch will be computed later in GLM_Init. ------------------------------------------------------- */ #ifdef GLM_MHD dxmin = grid[IDIR].dl_min; for (idim = 1; idim < DIMENSIONS; idim++){ /* Min cell length */ dxmin = MIN(dxmin, grid[idim].dl_min); } glm_ch = ini->cfl*dxmin/dtnext; #endif /* --------------------------------------------------------- 4. With STS, the ratio between advection (full) and parabolic time steps should not exceed ini->rmax_par. --------------------------------------------------------- */ #if (PARABOLIC_FLUX & SUPER_TIME_STEPPING) || (PARABOLIC_FLUX & RK_CHEBYSHEV) dt_par = ini->cfl_par/(2.0*Dts->inv_dtp); dtnext *= MIN(1.0, ini->rmax_par/(dt_adv/dt_par)); #endif /* ---------------------------------- 5. Compute Cooling time step ---------------------------------- */ #if COOLING != NO dtnext = MIN(dtnext, Dts->dt_cool); #endif /* -------------------------------------------------------------- 6. Allow time step to vary at most by a factor ini->cfl_max_var. Quit if dt gets too small, issue a warning if first_dt has been overestimated. -------------------------------------------------------------- */ dtnext = MIN(dtnext, ini->cfl_max_var*g_dt); if (dtnext < ini->first_dt*1.e-9){ print1 ("! dt is too small (%12.6e)!\n", dtnext); print1 ("! Cannot continue\n"); QUIT_PLUTO(1); } if (g_stepNumber <= 1 && (ini->first_dt > dtnext/ini->cfl)){ print1 ("! NextTimeStep: initial dt exceeds stability limit\n"); } /* -------------------------------------------- 7. Reset time step coefficients -------------------------------------------- */ DIM_LOOP(idim) Dts->cmax[idim] = 0.0; Dts->inv_dta = 0.0; Dts->inv_dtp = 0.0; Dts->dt_cool = 1.e38; return(dtnext); }
/* ********************************************************************* */ void EntropyOhmicHeating (const Data *d, Data_Arr UU, double dt, Grid *grid) /*! * Add Ohmic heating term to the conservative entropy equation when * RESISTIVITY is set to YES. * * \return This function has no return value. *********************************************************************** */ { #if PHYSICS == MHD && (RESISTIVITY == EXPLICIT) /* at the moment ... remove later ! */ int i,j,k,nv; double rho, rhog, gm1, vc[NVAR]; double Jc[3], eta[3], J2eta; double ***Jx1, *x1, *dx1; double ***Jx2, *x2, *dx2; double ***Jx3, *x3, *dx3; Data_Arr V; /* --------------------------------------------- 1. Set pointer shortcuts --------------------------------------------- */ V = d->Vc; Jx1 = d->J[IDIR]; x1 = grid[IDIR].x ; dx1 = grid[IDIR].dx; Jx2 = d->J[JDIR]; x2 = grid[JDIR].x ; dx2 = grid[JDIR].dx; Jx3 = d->J[KDIR]; x3 = grid[KDIR].x ; dx3 = grid[KDIR].dx; gm1 = g_gamma - 1.0; Jc[IDIR] = Jc[JDIR] = Jc[KDIR] = 0.0; /* --------------------------------------------- 2. Main spatial loop --------------------------------------------- */ DOM_LOOP(k,j,i){ /* --------------------------------- 2a. compute currents at cell center --------------------------------- */ #ifdef STAGGERED_MHD /* Staggered MHD */ #if COMPONENTS == 3 Jc[IDIR] = AVERAGE_YZ(Jx1,k-1,j-1,i); Jc[JDIR] = AVERAGE_XZ(Jx2,k-1,j,i-1); #endif Jc[KDIR] = AVERAGE_XY(Jx3,k,j-1,i-1); #else /* Cell-centered MHD */ #if COMPONENTS == 3 Jc[IDIR] = CDIFF_X2(V[BX3],k,j,i)/dx2[j] - CDIFF_X3(V[BX2],k,j,i)/dx3[k]; Jc[JDIR] = CDIFF_X3(V[BX1],k,j,i)/dx3[k] - CDIFF_X1(V[BX3],k,j,i)/dx1[i]; #endif Jc[KDIR] = CDIFF_X1(V[BX2],k,j,i)/dx1[i] - CDIFF_X2(V[BX1],k,j,i)/dx2[j]; #if GEOMETRY != CARTESIAN print1 ("! EntropyOhmicHeating: only CT supported in this geometry.\n") QUIT_PLUTO(1); #endif #endif /* ---------------------------------------- 2b. compute resistivity at cell center. ---------------------------------------- */ VAR_LOOP(nv) vc[nv] = V[nv][k][j][i]; rho = vc[RHO]; rhog = pow(rho,-gm1); Resistive_eta (vc, x1[i], x2[j], x3[k], Jc, eta); J2eta = 0.0; for (nv = 0; nv < 3; nv++) J2eta += Jc[nv]*Jc[nv]*eta[nv]; /* ---------------------------------------- 2c. Update conserved entropy ---------------------------------------- */ UU[k][j][i][ENTR] += dt*rhog*gm1*J2eta; }
/* ********************************************************************* */ void ParseCmdLineArgs (int argc, char *argv[], char *ini_file, Cmd_Line *cmd) /*! * * Parse command line options. Error messages will be output to * stdout. * * \param[in] argc argument count * \param[in] argv argument vector * \param[in] ini_file a pointer to a string containing the name of * the initialization file (default: "pluto.ini") * \param[out] cmd the command-line structure. * *********************************************************************** */ { int i,j; /* ----------------------------------------------- Set default values first ----------------------------------------------- */ cmd->restart = NO; cmd->h5restart = NO; cmd->maxsteps = 0 ; cmd->write = YES; cmd->makegrid = NO; cmd->jet = -1; /* -- means no direction -- */ cmd->show_dec = NO; cmd->xres = -1; /* -- means no grid resizing -- */ /* AYW -- 2012-06-19 09:21 JST maxtime default - 0 means no limit set*/ cmd->maxtime = 0 ; /* --AYW */ /* DM 10Aug15: initialise only */ cmd->init_only = NO; cmd->nproc[IDIR] = -1; /* means autodecomp will be used */ cmd->nproc[JDIR] = -1; cmd->nproc[KDIR] = -1; #ifdef PARALLEL cmd->parallel_dim[IDIR] = YES; /* by default, we parallelize */ cmd->parallel_dim[JDIR] = YES; /* all directions */ cmd->parallel_dim[KDIR] = YES; #endif /* -------------------------------------------------------------------- Parse Command Line Options. Note: Since at this time the output directory is now known and "pluto.log" has not been opened yet, we use printf to issue error messages. -------------------------------------------------------------------- */ for (i = 1; i < argc ; i++){ if (!strcmp(argv[i],"-dec")) { /* -- start reading integers at i+1 -- */ for (g_dir = 0; g_dir < DIMENSIONS; g_dir++){ if ((++i) >= argc){ if (prank == 0){ D_SELECT(printf ("! You must specify -dec n1\n"); , printf ("! You must specify -dec n1 n2\n"); , printf ("! You must specify -dec n1 n2 n3\n");) } QUIT_PLUTO(1); } cmd->nproc[g_dir] = atoi(argv[i]); if (cmd->nproc[g_dir] == 0){ if (prank == 0) { printf ("! Incorrect number of processor for g_dir = %d \n", g_dir); } QUIT_PLUTO(0); } }
/* ********************************************************************* */ void Init (double *us, double x1, double x2, double x3) /* * * * *********************************************************************** */ { double dr, vol, r; g_gamma = g_inputParam[GAMMA]; /* -------------------------------------------------- dr is the size of the initial energy deposition region: 2 ghost zones. -------------------------------------------------- */ #if DIMENSIONS == 1 dr = 2.0*(g_domEnd[IDIR]-g_domBeg[IDIR])/(double)NX1; #else dr = 3.5*(g_domEnd[IDIR]-g_domBeg[IDIR])/(double)NX1; #endif /* --------------------------------------- compute region volume --------------------------------------- */ #if (GEOMETRY == CARTESIAN) && (DIMENSIONS == 1) vol = 2.0*dr; #elif (GEOMETRY == CYLINDRICAL && DIMENSIONS == 1)|| \ (GEOMETRY == CARTESIAN && DIMENSIONS == 2) vol = CONST_PI*dr*dr; #elif (GEOMETRY == SPHERICAL && DIMENSIONS == 1)|| \ (GEOMETRY == CYLINDRICAL && DIMENSIONS == 2)|| \ (GEOMETRY == CARTESIAN && DIMENSIONS == 3) vol = 4.0/3.0*CONST_PI*dr*dr*dr; #else print1 ("! Init: geometrical configuration not allowed\n"); QUIT_PLUTO(1); #endif r = EXPAND(x1*x1, + x2*x2, +x3*x3); r = sqrt(r); us[RHO] = g_inputParam[DNST0]; us[VX1] = 0.0; us[VX2] = 0.0; us[VX3] = 0.0; if (r <= dr) { us[PRS] = (g_gamma - 1.0)*g_inputParam[ENRG0]/vol; }else{ us[PRS] = 1.e-5; } us[TRC] = 0.0; #if PHYSICS == MHD || PHYSICS == RMHD us[BX1] = 0.0; us[BX2] = 0.0; us[BX3] = 0.0; #ifdef STAGGERED_MHD us[AX] = 0.0; us[AY] = 0.0; us[AZ] = 0.0; #endif #endif }
/* ********************************************************************* */ int Setup (Input *input, Cmd_Line *cmd_line, char *ini_file) /*! * Open and parse the initialization file. * Assign values to the input structure. * * \param [out] input pointer to an Input structure * \param [in] cmd_line pointer to a Cmd_Line structure (useful, e.g., * to resize the domain using the \c -xres option) * \param [in] ini_file the name of the initialization file (default * is "pluto.ini") specified with the \c -i option. * *********************************************************************** */ { int idim, ip, ipos, itype, nlines; char *bound_opt[NOPT], str_var[512], *str; char *glabel[] = {"X1-grid", "X2-grid","X3-grid"}; char *bbeg_label[] = {"X1-beg", "X2-beg","X3-beg"}; char *bend_label[] = {"X1-end", "X2-end","X3-end"}; double dbl_var, rx; Output *output; FILE *fp; print1 ("> Reading %s (Setup) ...\n\n", ini_file); for (itype = 0; itype < NOPT; itype++) { bound_opt[itype] = "0000"; } /* --------------------------------------------------- available options are given as two set of names; This facilitates when updating the code and people are too lazy to read the manual ! --------------------------------------------------- */ bound_opt[OUTFLOW] = "outflow"; bound_opt[REFLECTIVE] = "reflective"; bound_opt[AXISYMMETRIC] = "axisymmetric"; bound_opt[EQTSYMMETRIC] = "eqtsymmetric"; bound_opt[PERIODIC] = "periodic"; bound_opt[SHEARING] = "shearingbox"; bound_opt[USERDEF] = "userdef"; input->log_freq = 1; /* -- default -- */ nlines = ParOpen (ini_file); /* ------------------------------------------------------------ [Grid] Section ------------------------------------------------------------ */ for (idim = 0; idim < 3; idim++){ input->npatch[idim] = atoi(ParGet(glabel[idim], 1)); input->npoint[idim] = 0; ipos = 1; for (ip = 1; ip <= input->npatch[idim]; ip++) { input->patch_left_node[idim][ip] = atof(ParGet(glabel[idim], ++ipos)); input->patch_npoint[idim][ip] = atoi(ParGet(glabel[idim], ++ipos)); input->npoint[idim] += input->patch_npoint[idim][ip]; input->grid_is_uniform[idim] = 0; strcpy (str_var, ParGet(glabel[idim], ++ipos)); /* printf ("%f %d %s\n",input->patch_left_node[idim][ip],input->patch_npoint[idim][ip],str_var); */ if (strcmp(str_var,"u") == 0 || strcmp(str_var,"uniform") == 0) { input->patch_type[idim][ip] = UNIFORM_GRID; if (input->npatch[idim] == 1) input->grid_is_uniform[idim] = 1; }else if (strcmp(str_var,"s") == 0 || strcmp(str_var,"strecthed") == 0) { input->patch_type[idim][ip] = STRETCHED_GRID; }else if (strcmp(str_var,"l+") == 0){ input->patch_type[idim][ip] = LOGARITHMIC_INC_GRID; }else if (strcmp(str_var,"l-") == 0){ input->patch_type[idim][ip] = LOGARITHMIC_DEC_GRID; }else{ print("\nYou must specify either 'u', 's', 'l+' or 'l-' as grid-type in %s\n", ini_file); QUIT_PLUTO(1); } } input->patch_left_node[idim][ip] = atof(ParGet(glabel[idim], ++ipos)); if ( (ipos+1) != (input->npatch[idim]*3 + 3)) { print (" ! Domain #%d setup is not properly defined \n", idim); QUIT_PLUTO(1); } if (idim >= DIMENSIONS && input->npoint[idim] != 1) { print ("! %d point(s) on dim. %d is NOT valid, resetting to 1\n", input->npoint[idim],idim+1); input->npoint[idim] = 1; input->npatch[idim] = 1; input->patch_npoint[idim][1] = 1; } } /* ------------------------------------------------------------ Change the resolution if cmd_line->xres has been given ------------------------------------------------------------ */ if (cmd_line->xres > 1) { rx = (double)cmd_line->xres/(double)input->patch_npoint[IDIR][1]; for (idim = 0; idim < DIMENSIONS; idim++){ if (input->npatch[idim] > 1){ print ("! -xres option works on uniform, single patch grid\n"); QUIT_PLUTO(1); } dbl_var = (double)input->patch_npoint[idim][1]; input->patch_npoint[idim][1] = (int)(dbl_var*rx); dbl_var = (double)input->npoint[idim]; input->npoint[idim] = (int)(dbl_var*rx); } } /* ------------------------------------------------------------ [Time] Section ------------------------------------------------------------ */ input->cfl = atof(ParGet("CFL", 1)); if (ParQuery ("CFL_par")) input->cfl_par = atof(ParGet("CFL_par", 1)); else input->cfl_par = 0.8/(double)DIMENSIONS; if (ParQuery ("rmax_par")) input->rmax_par = atof(ParGet("rmax_par", 1)); else input->rmax_par = 100.0; input->cfl_max_var = atof(ParGet("CFL_max_var", 1)); input->tstop = atof(ParGet("tstop", 1)); input->first_dt = atof(ParGet("first_dt", 1)); /* ------------------------------------------------------------ [Solver] Section ------------------------------------------------------------ */ sprintf (input->solv_type,"%s",ParGet("Solver",1)); /* ------------------------------------------------------------ [Boundary] Section ------------------------------------------------------------ */ for (idim = 0; idim < 3; idim++){ str = ParGet(bbeg_label[idim], 1); COMPARE (str, bound_opt[itype], itype); if (itype == NOPT) { print ("! Setup: don't know how to put left boundary '%s' \n", str); QUIT_PLUTO(1); } input->lft_bound_side[idim] = itype; } for (idim = 0; idim < 3; idim++){ str = ParGet(bend_label[idim], 1); COMPARE (str, bound_opt[itype], itype); if (itype == NOPT) { print ("! Setup: don't know how to put left boundary '%s' \n", str); QUIT_PLUTO(1); } input->rgt_bound_side[idim] = itype; } /* ------------------------------------------------------------ [Output] Section ------------------------------------------------------------ */ input->user_var = atoi(ParGet("uservar", 1)); for (ip = 0; ip < input->user_var; ip++){ if ( (str = ParGet("uservar", 2 + ip)) != NULL){ sprintf (input->user_var_name[ip], "%s", str); }else{ print ("! Setup: missing name after user var name '%s'\n", input->user_var_name[ip-1]); QUIT_PLUTO(1); } } /* ---- dbl output ---- */ ipos = 0; output = input->output + (ipos++); output->type = DBL_OUTPUT; GetOutputFrequency(output, "dbl"); sprintf (output->mode,"%s",ParGet("dbl",3)); #ifdef USE_ASYNC_IO if ( strcmp(output->mode,"single_file") && strcmp(output->mode,"single_file_async") && strcmp(output->mode,"multiple_files")){ print1 ( "! Setup: expecting 'single_file', 'single_file_async' or 'multiple_files' in dbl output\n"); QUIT_PLUTO(1); } #else if ( strcmp(output->mode,"single_file") && strcmp(output->mode,"multiple_files")){ print1 ( "! Setup: expecting 'single_file' or 'multiple_files' in dbl output\n"); QUIT_PLUTO(1); } #endif /* ---- flt output ---- */ output = input->output + (ipos++); output->type = FLT_OUTPUT; GetOutputFrequency(output, "flt"); sprintf (output->mode,"%s",ParGet("flt",3)); #ifdef USE_ASYNC_IO if ( strcmp(output->mode,"single_file") && strcmp(output->mode,"single_file_async") && strcmp(output->mode,"multiple_files")){ print1 ( "! Setup: expecting 'single_file', 'single_file_async' or 'multiple_files' in flt output\n"); QUIT_PLUTO(1); } #else if ( strcmp(output->mode,"single_file") && strcmp(output->mode,"multiple_files")){ print1 ( "! Setup: expecting 'single_file' or 'multiple_files' in flt output\n"); QUIT_PLUTO(1); } #endif /* -- hdf5 output -- */ if (ParQuery("dbl.h5")){ output = input->output + (ipos++); output->type = DBL_H5_OUTPUT; GetOutputFrequency(output, "dbl.h5"); } if (ParQuery("flt.h5")){ output = input->output + (ipos++); output->type = FLT_H5_OUTPUT; GetOutputFrequency(output, "flt.h5"); } /* -- vtk output -- */ if (ParQuery ("vtk")){ output = input->output + (ipos++); output->type = VTK_OUTPUT; GetOutputFrequency(output, "vtk"); if (ParGet("vtk",3) == NULL){ print1 (" ! Setup: extra field missing in vtk output\n"); QUIT_PLUTO(1); } sprintf (output->mode,"%s",ParGet("vtk",3)); if ( strcmp(output->mode,"single_file") && strcmp(output->mode,"multiple_files")){ print1 (" ! Setup: expecting 'single_file' or 'multiple_files' in\n"); print1 (" vtk output\n"); QUIT_PLUTO(1); } } /* -- tab output -- */ if (ParQuery ("tab")){ output = input->output + (ipos++); output->type = TAB_OUTPUT; GetOutputFrequency(output, "tab"); } /* -- ppm output -- */ if (ParQuery ("ppm")){ output = input->output + (ipos++); output->type = PPM_OUTPUT; GetOutputFrequency(output, "ppm"); } /* -- png output -- */ if (ParQuery ("png")){ output = input->output + (ipos++); output->type = PNG_OUTPUT; GetOutputFrequency(output, "png"); } /* -- log frequency -- */ input->log_freq = atoi(ParGet("log", 1)); input->log_freq = MAX(input->log_freq, 1); /* -- set default for remaining output type -- */ while (ipos < MAX_OUTPUT_TYPES){ output = input->output + ipos; output->type = -1; output->dt = -1.0; output->dn = -1; output->dclock = -1.0; ipos++; } /* -- analysis -- */ if (ParQuery ("analysis")){ input->anl_dt = atof(ParGet("analysis", 1)); input->anl_dn = atoi(ParGet("analysis", 2)); }else{ input->anl_dt = -1.0; /* -- defaults -- */ input->anl_dn = -1; } /* ------------------------------------------------------------ [Parameters] Section ------------------------------------------------------------ */ fp = fopen(ini_file,"r"); /* -- find position at "[Parameters" -- */ for (ipos = 0; ipos <= nlines; ipos++){ fgets(str_var, 512, fp); if (strlen(str_var) > 0) { str = strtok (str_var,"]"); if (strcmp(str,"[Parameters") == 0) break; } } fgets(str_var, 512, fp); for (ip = 0; ip < USER_DEF_PARAMETERS; ip++){ fscanf (fp,"%s %lf\n", str_var, &dbl_var); input->aux[ip] = dbl_var; } fclose(fp); /* -------------------------------------------- print some output -------------------------------------------- */ print1 (" CFL : %4.2f\n",input->cfl); print1 (" Solver: %s\n",input->solv_type); for (ip = 0; ip < USER_DEF_PARAMETERS; ip++){ print1 (" User def par #%d = %10.4e\n",ip,input->aux[ip]); } for (idim = 0; idim < DIMENSIONS; idim++){ print1 (" X%d boundary: ",idim+1); print1 ("[beg ... end] %s ... %s\n", bound_opt[input->lft_bound_side[idim]], bound_opt[input->rgt_bound_side[idim]]); } return(0); }
/* ********************************************************************* */ void Init (double *v, double x, double y, double z) /* * *********************************************************************** */ { static int first_call = 1; double Bz0, rnd, dvy, cs, H; double Lx, Ly, Lz; double kx, ky, kz; #ifndef SHEARINGBOX print1 ("! ShearingBox module has not been included.\n"); print1 ("! Cannot continue.\n"); QUIT_PLUTO(1); #endif /* -- compute domain sizes -- */ Lx = g_domEnd[IDIR] - g_domBeg[IDIR]; kx = 2.0*CONST_PI/Lx; Ly = g_domEnd[JDIR] - g_domBeg[JDIR]; ky = 2.0*CONST_PI/Ly; Lz = g_domEnd[KDIR] - g_domBeg[KDIR]; kz = 2.0*CONST_PI/Lz; /* -- get sound speed and pressure scale height -- */ cs = g_inputParam[CSOUND]; /* sound speed */ //H = sqrt(2.0)*cs/sb_Omega; /* pressure scale height */ H = sqrt(2.0)*cs/SB_OMEGA; /* pressure scale height */ /* -- seed random numbers differently for each processor -- */ if (first_call == 1){ srand(time(NULL) + prank); first_call = 0; } rnd = (double)(rand())/((double)RAND_MAX + 1.0); /* -- set velocity perturbation [1]: random noise -- */ /* dvy = 0.01*cs*rnd; */ /* -- set velocity perturbation [2], two harmonics for each direction -- */ dvy = sin(kx*x + 0.20) + sin(2.0*kx*x - 0.37); dvy *= sin(ky*y + 0.13) + sin(2.0*ky*y + 0.04); dvy *= sin(kz*z + 0.56) + sin(2.0*kz*z + 0.62); dvy *= 0.01*cs/8.0; /* -- in 2D we don't use any perturbation -- */ #if DIMENSIONS == 2 dvy = 0.0; #endif /* -- set initial condition -- */ #if STRATIFICATION == YES v[RHO] = exp(-z*z/(H*H)); #else v[RHO] = 1.0; #endif v[VX1] = 0.0; //v[VX2] = -sb_q*sb_Omega*x + dvy; v[VX2] = -SB_Q*SB_OMEGA*x + dvy; v[VX3] = 0.0; #if EOS == IDEAL v[PRS] = cs*cs*v[RHO]; #elif EOS == ISOTHERMAL g_isoSoundSpeed = cs; #endif v[TRC] = 0.0; /* ---------------------------------------------------------------- The magnetic field amplitude is set by the parameter beta = 2p/B^2 = 2*rho*c^2/b0^2 --> b0 = c*sqrt(2/beta) where it is assumed that rho = 1 (midplane). ---------------------------------------------------------------- */ #if PHYSICS == MHD Bz0 = cs*sqrt(2.0/g_inputParam[BETA]); #if NET_FLUX == YES /* -- Net flux, constant vertical field Bz = B0 -- */ v[BX1] = 0.0; v[BX2] = 0.0; v[BX3] = Bz0; v[AX1] = 0.0; v[AX2] = Bz0*x; v[AX3] = 0.0; #else /* -- Zero net flux, Bz = B0*sin(2*pi*x/Lx) -- */ v[BX1] = 0.0; v[BX2] = 0.0; v[BX3] = Bz0*sin(kx*x); v[AX1] = 0.0; v[AX2] = -Bz0*cos(kx*x)/kx; v[AX3] = 0.0; #endif #if DIMENSIONS == 2 /* 2D Case only for testing */ v[BX1] = Bz0*sin(ky*y); v[BX2] = 0.0; v[BX3] = 0.0; v[AX1] = 0.0; v[AX2] = 0.0; v[AX3] = -Bz0*cos(ky*y)/ky; #endif #endif }
/* **************************************************************** */ void SoundSpeed2 (double **u, double *cs2, double *h, int beg, int end, int pos, Grid *grid) /* * * Define the square of the sound speed for different EOS * ****************************************************************** */ { int i; #if EOS == IDEAL for (i = beg; i <= end; i++) cs2[i] = g_gamma*u[i][PRS]/u[i][RHO]; #elif EOS == ISOTHERMAL { int j,k; /* -- used as multidimensional indices -- */ double *x1, *x2, *x3; x1 = grid[IDIR].x; x2 = grid[JDIR].x; x3 = grid[KDIR].x; i = g_i; j = g_j; k = g_k; if (g_dir == IDIR) { double R; x1 = (pos == FACE_CENTER ? grid[IDIR].xr : grid[IDIR].x); for (i = beg; i <= end; i++){ #if GEOMETRY == POLAR R = x1[i]; #elif GEOMETRY == SPHERICAL R = x1[i]*sin(x2[j]); #endif cs2[i] = g_isoSoundSpeed*g_isoSoundSpeed/R; } }else if (g_dir == JDIR){ double R; x2 = (pos == FACE_CENTER ? grid[JDIR].xr : grid[JDIR].x); for (j = beg; j <= end; j++) { #if GEOMETRY == POLAR R = x1[i]; #elif GEOMETRY == SPHERICAL R = x1[i]*sin(x2[j]); #endif cs2[j] = g_isoSoundSpeed*g_isoSoundSpeed/R; } }else if (g_dir == KDIR){ double R; x3 = (pos == FACE_CENTER ? grid[KDIR].xr : grid[KDIR].x); for (k = beg; k <= end; k++){ #if GEOMETRY == POLAR R = x1[i]; #elif GEOMETRY == SPHERICAL R = x1[i]*sin(x2[j]); #endif cs2[k] = g_isoSoundSpeed*g_isoSoundSpeed/R; } } } #else print ("! SoundSpeed2: not defined for this EoS\n"); QUIT_PLUTO(1); #endif }
/* ********************************************************************* */ void FARGO_ShiftSolution(Data_Arr U, Data_Arr Us, Grid *grid) /*! * Shifts conserved variables in the orbital direction. * * \param [in,out] U a 3D array of conserved, zone-centered values * \param [in,out] Us a 3D array of staggered magnetic fields * \param [in] grid pointer to Grid structure; * * \return This function has no return value. * \todo Optimization: avoid using too many if statement like * on nproc_s > 1 or == 1 *********************************************************************** */ #if GEOMETRY != SPHERICAL #define s j /* -- orbital direction index -- */ #else #define s k #endif { int i, j, k, nv, ngh, nvar_c; int sm, m; static int mmax, mmin; /* Used to determine buffer size from time to time */ int nbuf_lower, nbuf_upper, nproc_s = 1; double *dx, *dy, *dz, dphi; double *x, *xp, *y, *yp, *z, *zp, w, Lphi, dL, eps; double **wA, ***Uc[NVAR]; static double *q00, *flux00; double *q, *flux; static double ***Ez, ***Ex; #ifdef PARALLEL double ***upper_buf[NVAR]; /* upper layer buffer */ double ***lower_buf[NVAR]; /* lower layer buffer */ RBox lbox, ubox; /* lower and upper box structures */ #endif #if GEOMETRY == CYLINDRICAL print1 ("! FARGO cannot be used in this geometry\n"); QUIT_PLUTO(1); #endif /* ----------------------------------------------------------------- Allocate memory for static arrays. In parallel mode, one-dimensional arrays must be large enough to contain data values coming from neighbour processors. For this reason we augment them with an extra zone buffer containing MAX_BUF_SIZE cells on both sides. In this way, the array indices for q can be negative. <..MAX_BUF_SIZE..>|---|---| ... |---|----|<..MAX_BUF_SIZE..> 0 1 NX2-1 ----------------------------------------------------------------- */ if (q00 == NULL){ #if GEOMETRY == SPHERICAL && DIMENSIONS == 3 q00 = ARRAY_1D(NX3_TOT + 2*MAX_BUF_SIZE,double); flux00 = ARRAY_1D(NX3_TOT + 2*MAX_BUF_SIZE,double); #else q00 = ARRAY_1D(NX2_TOT + 2*MAX_BUF_SIZE,double); flux00 = ARRAY_1D(NX2_TOT + 2*MAX_BUF_SIZE,double); #endif #ifdef STAGGERED_MHD Ez = ARRAY_3D(NX3_TOT, NX2_TOT, NX1_TOT, double); Ex = ARRAY_3D(NX3_TOT, NX2_TOT, NX1_TOT, double); #endif /* -------------------------------------------------------------- The value of mmax and mmin is kept from call to call because they are used to determine the size of the send/receive buffers each time step. At the very beginning, we initialize them to the largest possible value. -------------------------------------------------------------- */ mmax = MIN(NS, MAX_BUF_SIZE-1) - grid[SDIR].nghost - 2; mmin = mmax; }
/* ********************************************************************* */ void FD_Flux (const State_1D *state, int beg, int end, double *cmax, Grid *grid) /*! * Compute interface flux by suitable high-order finite * difference non-oscillatory interpolants. * * \param [in] state pointer to State_1D structure * \param [in] beg initial index of computation * \param [in] end final index of computation * \param [out] cmax array of maximum characteristic speeds * \param [in] grid pointer to an array of Grid structures * * \b Reference * - "High-order conservative finite difference GLM-MHD schemes * for cell-centered MHD"\n * Mignone et al., JCP (2010), 229, 5896 * * \return This function has no return value. * ****************************************************************** */ { int nv, i, j, k, np_tot; int Kmax, S; double **v, **flux, *press; double fs, fp, fm, flx[NFLX], dx; static double *Fp, *Fm, **F, *a2, **u; static double **Uave, **Vave, **lp; double **L, **R; static double *psim, *Bm; double (*REC)(double *, double, int); /* ------------------------------------------------------------------ Check compatibility with other modules ------------------------------------------------------------------ */ #if TIME_STEPPING != RK3 print1 ("! Finite Difference schemes work with RK3 only \n"); QUIT_PLUTO(1); #endif #if GEOMETRY != CARTESIAN print1 ("! Finite Difference schemes work in Cartesian coordinates only\n"); QUIT_PLUTO(1); #endif #if PHYSICS == RMHD || PHYSICS == RHD print1 ("! Finite difference schemes work only for HD od MHD modules\n"); QUIT_PLUTO(1); #endif #ifdef STAGGERED_MHD print1 ("! Finite difference schemes work only with cell-centered schemes\n"); QUIT_PLUTO(1); #endif /* -------------------------------------------------------------- - Define pointer to reconstruction function - Determine the stencil for interpolation For a given S, the stencil is: i-S <= i <= i+S+1 -------------------------------------------------------------- */ dx = grid[g_dir].dx[beg]; #if INTERPOLATION == WENO3_FD REC = WENO3_Reconstruct; S = 1; #elif INTERPOLATION == LIMO3_FD REC = LIMO3_Reconstruct; S = 1; #elif INTERPOLATION == WENOZ_FD REC = WENOZ_Reconstruct; S = 2; #elif INTERPOLATION == MP5_FD REC = MP5_Reconstruct; S = 2; #endif if (F == NULL){ Fp = ARRAY_1D(NMAX_POINT, double); Fm = ARRAY_1D(NMAX_POINT, double); F = ARRAY_2D(NMAX_POINT, NVAR, double); lp = ARRAY_2D(NMAX_POINT, NVAR, double); psim = ARRAY_1D(NMAX_POINT, double); Bm = ARRAY_1D(NMAX_POINT, double); a2 = ARRAY_1D(NMAX_POINT, double); Uave = ARRAY_2D(NMAX_POINT, NVAR, double); Vave = ARRAY_2D(NMAX_POINT, NVAR, double); u = ARRAY_2D(NMAX_POINT, NVAR, double); }
/* ********************************************************************* */ double GaussQuadrature(double (*func)(double), double xb, double xe, int nstep, int order) /*! * Perform numerical quadrature of the function f(x) between * the lower bound xb and upper bound xe by subdividing the interval * into 'nstep' steps. * A 3 or 5-point Gaussian quadrature rule is used depending on the * input variable order (=3 or =5) * * \param [in] *func a pointer to the function func(x) (returning double) * to be integrated * \param [in] xb the lower interval bound * \param [in] xe the upper interval bound * \param [in] nstep the number of sub-intervals into which the * original interval [xb,xe] has to be divided * \param [in] order the number of Gaussian points (only 3 or 5) * *********************************************************************** */ { int i, n; double w[8], z[8], x; double I, Isub; double xb0, xe0, dx; if (order == 3){ double s3 = sqrt(3.0/5.0); z[0] = -s3; w[0] = 5.0/9.0; z[1] = 0.0; w[1] = 8.0/9.0; z[2] = s3; w[2] = 5.0/9.0; }else if (order == 5){ double s1, s7; s1 = sqrt(10.0/7.0); s7 = sqrt(70.0); z[0] = -1.0/3.0*sqrt(5.0 - 2.0*s1); w[0] = (322.0 + 13.0*s7)/900.0; z[1] = -z[0]; w[1] = w[0]; z[2] = -1.0/3.0*sqrt(5.0 + 2.0*s1); w[2] = (322.0 - 13.0*s7)/900.0; z[3] = -z[2]; w[3] = w[2]; z[4] = 0.0; w[4] = 128.0/225.0; }else{ print ("! GaussQuadrature: order must be either 3 or 5\n"); QUIT_PLUTO(1); } if (nstep <= 0){ print ("! GaussQuadrature: nstep must be > 0\n"); QUIT_PLUTO(1); } xb0 = xb; xe0 = xe; /* save original interval endpoints */ dx = (xe - xb)/(double)nstep; /* sub-interval length */ I = 0.0; for (i = 0; i < nstep; i++){ xb = xb0 + i*dx; xe = xb + dx; Isub = 0.0; /* intgrate sub-interval */ for (n = 0; n < order; n++){ x = 0.5*(xe - xb)*z[n] + (xe + xb)*0.5; Isub += w[n]*func(x); } Isub *= 0.5*(xe - xb); I += Isub; } return I; }
/* ********************************************************************* */ void WriteData (const Data *d, Output *output, Grid *grid) /*! * Write data to disk using any of the available formats. * * \param [in] d pointer to PLUTO Data structre * \param [in] output the output structure corresponding to a given * format * \param [in] grid pointer to an array of Grid structures *********************************************************************** */ { int i, j, k, nv; int single_file; size_t dsize; char filename[512], sline[512]; static int last_computed_var = -1; double units[MAX_OUTPUT_VARS]; float ***Vpt3; void *Vpt; FILE *fout, *fbin; time_t tbeg, tend; long long offset; /* ----------------------------------------------------------- Increment the file number and initialize units ----------------------------------------------------------- */ output->nfile++; print1 ("> Writing file #%d (%s) to disk...", output->nfile, output->ext); #ifdef PARALLEL MPI_Barrier (MPI_COMM_WORLD); if (prank == 0) time(&tbeg); #endif for (nv = 0; nv < MAX_OUTPUT_VARS; nv++) units[nv] = 1.0; if (output->cgs) GetCGSUnits(units); /* -------------------------------------------------------- Get user var if necessary -------------------------------------------------------- */ if (last_computed_var != g_stepNumber && d->Vuser != NULL) { ComputeUserVar (d, grid); last_computed_var = g_stepNumber; } /* -------------------------------------------------------- Select the output type -------------------------------------------------------- */ if (output->type == DBL_OUTPUT) { /* ------------------------------------------------------------------- */ /*! - \b DBL output: Double-precision data files can be written using single or multiple file mode. - for single file, serial: we open the file just once before the main variable loop, dump variables and then close. - for single file, parallel the distributed array descriptor sz is different for cell-centered or staggered data type and we thus have to open and close the file after each variable has been dumped. - when writing multiple files we open, write to and close the file one each loop cycle. \note In all cases, the pointer to the data array that has to be written must be cast into (void *) and the starting index of the array must be zero. */ /* ------------------------------------------------------------------- */ int sz; single_file = strcmp(output->mode,"single_file") == 0; dsize = sizeof(double); if (single_file){ /* -- single output file -- */ sprintf (filename, "%s/data.%04d.%s", output->dir,output->nfile, output->ext); offset = 0; #ifndef PARALLEL fbin = OpenBinaryFile (filename, 0, "w"); #endif for (nv = 0; nv < output->nvar; nv++) { if (!output->dump_var[nv]) continue; if (output->stag_var[nv] == -1) { /* -- cell-centered data -- */ sz = SZ; Vpt = (void *)output->V[nv][0][0]; } else if (output->stag_var[nv] == 0) { /* -- x-staggered data -- */ sz = SZ_stagx; Vpt = (void *)(output->V[nv][0][0]-1); } else if (output->stag_var[nv] == 1) { /* -- y-staggered data -- */ sz = SZ_stagy; Vpt = (void *)output->V[nv][0][-1]; } else if (output->stag_var[nv] == 2) { /* -- z-staggered data -- */ sz = SZ_stagz; Vpt = (void *)output->V[nv][-1][0]; } #ifdef PARALLEL fbin = OpenBinaryFile (filename, sz, "w"); AL_Set_offset(sz, offset); #endif WriteBinaryArray (Vpt, dsize, sz, fbin, output->stag_var[nv]); #ifdef PARALLEL offset = AL_Get_offset(sz); CloseBinaryFile(fbin, sz); #endif } #ifndef PARALLEL CloseBinaryFile(fbin, sz); #endif }else{ /* -- multiple files -- */ for (nv = 0; nv < output->nvar; nv++) { if (!output->dump_var[nv]) continue; sprintf (filename, "%s/%s.%04d.%s", output->dir, output->var_name[nv], output->nfile, output->ext); if (output->stag_var[nv] == -1) { /* -- cell-centered data -- */ sz = SZ; Vpt = (void *)output->V[nv][0][0]; } else if (output->stag_var[nv] == 0) { /* -- x-staggered data -- */ sz = SZ_stagx; Vpt = (void *)(output->V[nv][0][0]-1); } else if (output->stag_var[nv] == 1) { /* -- y-staggered data -- */ sz = SZ_stagy; Vpt = (void *)output->V[nv][0][-1]; } else if (output->stag_var[nv] == 2) { /* -- z-staggered data -- */ sz = SZ_stagz; Vpt = (void *)output->V[nv][-1][0]; } fbin = OpenBinaryFile (filename, sz, "w"); WriteBinaryArray (Vpt, dsize, sz, fbin, output->stag_var[nv]); CloseBinaryFile (fbin, sz); } } } else if (output->type == FLT_OUTPUT) { /* ---------------------------------------------------------- FLT output for cell-centered data ---------------------------------------------------------- */ single_file = strcmp(output->mode,"single_file") == 0; if (single_file){ /* -- single output file -- */ sprintf (filename, "%s/data.%04d.%s", output->dir, output->nfile, output->ext); fbin = OpenBinaryFile (filename, SZ_float, "w"); for (nv = 0; nv < output->nvar; nv++) { if (!output->dump_var[nv]) continue; /* Vpt = (void *)(Convert_dbl2flt(output->V[nv],0))[0][0]; */ Vpt3 = Convert_dbl2flt(output->V[nv], units[nv],0); Vpt = (void *)Vpt3[0][0]; WriteBinaryArray (Vpt, sizeof(float), SZ_float, fbin, output->stag_var[nv]); } CloseBinaryFile(fbin, SZ_float); /* BOV_Header(output, filename); */ }else{ /* -- multiple files -- */ for (nv = 0; nv < output->nvar; nv++) { if (!output->dump_var[nv]) continue; sprintf (filename, "%s/%s.%04d.%s", output->dir, output->var_name[nv], output->nfile, output->ext); fbin = OpenBinaryFile (filename, SZ_float, "w"); /* Vpt = (void *)(Convert_dbl2flt(output->V[nv],0))[0][0]; */ Vpt3 = Convert_dbl2flt(output->V[nv], units[nv],0); Vpt = (void *)Vpt3[0][0]; WriteBinaryArray (Vpt, sizeof(float), SZ_float, fbin, output->stag_var[nv]); CloseBinaryFile (fbin, SZ_float); } } }else if (output->type == DBL_H5_OUTPUT || output->type == FLT_H5_OUTPUT){ /* ------------------------------------------------------ HDF5 (static grid) output (single/double precision) ------------------------------------------------------ */ #ifdef USE_HDF5 single_file = YES; WriteHDF5 (output, grid); #else print1 ("! WriteData: HDF5 library not available\n"); return; #endif }else if (output->type == VTK_OUTPUT) { /* ------------------------------------------------------------------- */ /*! - \b VTK output: in order to enable parallel writing, files must be closed and opened again for scalars, since the distributed array descriptors used by ArrayLib (Float_Vect) and (float) are different. This is done using the AL_Get_offset() and AL_Set_offset() functions. */ /* ------------------------------------------------------------------- */ single_file = strcmp(output->mode,"single_file") == 0; sprintf (filename, "%s/data.%04d.%s", output->dir, output->nfile, output->ext); if (single_file){ /* -- single output file -- */ fbin = OpenBinaryFile(filename, SZ_Float_Vect, "w"); WriteVTK_Header(fbin, grid); for (nv = 0; nv < output->nvar; nv++) { /* -- write vectors -- */ if (output->dump_var[nv] != VTK_VECTOR) continue; WriteVTK_Vector (fbin, output->V + nv, units[nv], output->var_name[nv], grid); } #ifdef PARALLEL offset = AL_Get_offset(SZ_Float_Vect); CloseBinaryFile(fbin, SZ_Float_Vect); fbin = OpenBinaryFile(filename, SZ_float, "w"); AL_Set_offset(SZ_float, offset); #endif for (nv = 0; nv < output->nvar; nv++) { /* -- write scalars -- */ if (output->dump_var[nv] != YES) continue; WriteVTK_Scalar (fbin, output->V[nv], units[nv], output->var_name[nv], grid); } CloseBinaryFile(fbin, SZ_float); }else{ /* -- multiple output files -- */ for (nv = 0; nv < output->nvar; nv++) { /* -- write vectors -- */ if (output->dump_var[nv] != VTK_VECTOR) continue; if (strcmp(output->var_name[nv],"vx1") == 0) { sprintf (filename, "%s/vfield.%04d.%s", output->dir, output->nfile, output->ext); }else if (strcmp(output->var_name[nv],"bx1") == 0) { sprintf (filename, "%s/bfield.%04d.%s", output->dir, output->nfile, output->ext); }else{ print1 ("! WriteData: unknown vector type in VTK output\n"); QUIT_PLUTO(1); } fbin = OpenBinaryFile(filename, SZ_Float_Vect, "w"); WriteVTK_Header(fbin, grid); WriteVTK_Vector(fbin, output->V + nv, units[nv], output->var_name[nv], grid); CloseBinaryFile(fbin, SZ_Float_Vect); } for (nv = 0; nv < output->nvar; nv++) { /* -- write scalars -- */ if (output->dump_var[nv] != YES) continue; sprintf (filename, "%s/%s.%04d.%s", output->dir, output->var_name[nv], output->nfile, output->ext); fbin = OpenBinaryFile(filename, SZ_Float_Vect, "w"); WriteVTK_Header(fbin, grid); #ifdef PARALLEL offset = AL_Get_offset(SZ_Float_Vect); CloseBinaryFile(fbin, SZ_Float_Vect); fbin = OpenBinaryFile(filename, SZ_float, "w"); AL_Set_offset(SZ_float, offset); #endif WriteVTK_Scalar(fbin, output->V[nv], units[nv], output->var_name[nv], grid); CloseBinaryFile (fbin, SZ_float); } } }else if (output->type == TAB_OUTPUT) { /* ------------------------------------------------------ Tabulated (ASCII) output ------------------------------------------------------ */ single_file = YES; sprintf (filename,"%s/data.%04d.%s", output->dir, output->nfile, output->ext); WriteTabArray (output, filename, grid); }else if (output->type == PPM_OUTPUT) { /* ------------------------------------------------------ PPM output ------------------------------------------------------ */ single_file = NO; for (nv = 0; nv < output->nvar; nv++) { if (!output->dump_var[nv]) continue; sprintf (filename, "%s/%s.%04d.%s", output->dir, output->var_name[nv], output->nfile, output->ext); WritePPM (output->V[nv], output->var_name[nv], filename, grid); } }else if (output->type == PNG_OUTPUT) { /* ------------------------------------------------------ PNG output ------------------------------------------------------ */ #ifdef USE_PNG single_file = NO; for (nv = 0; nv < output->nvar; nv++) { if (!output->dump_var[nv]) continue; sprintf (filename, "%s/%s.%04d.%s", output->dir, output->var_name[nv], output->nfile, output->ext); WritePNG (output->V[nv], output->var_name[nv], filename, grid); } #else print1 ("! PNG library not available\n"); return; #endif } /* ------------------------------------------------------------- Update corresponding ".out" file ------------------------------------------------------------- */ sprintf (filename,"%s/%s.out",output->dir, output->ext); if (prank == 0) { if (output->nfile == 0) { fout = fopen (filename, "w"); }else { fout = fopen (filename, "r+"); for (nv = 0; nv < output->nfile; nv++) fgets (sline, 512, fout); fseek (fout, ftell(fout), SEEK_SET); } /* -- write a multi-column file -- */ fprintf (fout, "%d %12.6e %12.6e %ld ", output->nfile, g_time, g_dt, g_stepNumber); if (single_file) fprintf (fout,"single_file "); else fprintf (fout,"multiple_files "); if (IsLittleEndian()) fprintf (fout, "little "); else fprintf (fout, "big "); for (nv = 0; nv < output->nvar; nv++) { if (output->dump_var[nv]) fprintf (fout, "%s ", output->var_name[nv]); } fprintf (fout,"\n"); fclose (fout); } #ifdef PARALLEL MPI_Barrier (MPI_COMM_WORLD); if (prank == 0){ time(&tend); print1 (" [%5.2f sec]",difftime(tend,tbeg)); } #endif print1 ("\n"); }
/* ********************************************************************* */ double MeanMolecularWeight(double *v) /*! * * Return the mean molecular weight. * * \param [in] v array of primitive variables (including ions) * *********************************************************************** */ { int nv; double mu; #if COOLING == NO mu = (CONST_AH + FRAC_He*CONST_AHe + FRAC_Z*CONST_AZ) / (2.0 + FRAC_He + FRAC_Z*(1.0 + CONST_AZ*0.5)); #elif COOLING == TABULATED mu = (CONST_AH + FRAC_He*CONST_AHe + FRAC_Z*CONST_AZ) / (2.0 + FRAC_He + FRAC_Z*(1.0 + CONST_AZ*0.5)); #elif COOLING == SNEq mu = (CONST_AH + FRAC_He*CONST_AHe + FRAC_Z*CONST_AZ) / (2.0 + FRAC_He + 2.0*FRAC_Z - v[X_HI]); /* return ( (CONST_AH + frac_He*CONST_AHe + frac_Z*CONST_AZ) / (2.0 + frac_He + 2.0*frac_Z - v[X_HI])); */ #elif COOLING == H2_COOL double munum, muden; NIONS_LOOP(nv){ v[nv] = MAX(v[nv], 0.0); v[nv] = MIN(v[nv], 1.0); } v[X_H2] = MIN(v[X_H2], 0.5); double fn = v[X_HI]; double gn = v[X_H2]; double hn = v[X_HII]; mu = (CONST_AH + CONST_AHe*FRAC_He + CONST_AZ*FRAC_Z) / (fn + gn + 2*hn + FRAC_He + FRAC_Z + 0.5*CONST_AZ*FRAC_Z); /* double N_H = (H_MASS_FRAC/CONST_AH); double N_He = (He_MASS_FRAC/CONST_AHe); double N_Z = ((1.0 - H_MASS_FRAC - He_MASS_FRAC)/CONST_AZ); double fracHe = N_He/N_H; double fracZ = N_Z/N_H; double fn = v[X_HI]; double gn = v[X_H2]; double hn = v[X_HII]; munum = 1.0 + CONST_AHe*(fracHe) + CONST_AZ*(fracZ); muden = fn + gn + 2*hn + fracHe + fracZ + 0.5*CONST_AZ*(fracZ); return munum/muden; */ #elif COOLING == MINEq double mmw1, mmw2; int i, j; mmw1 = mmw2 = 0.0; for (i = 0; i < NIONS; i++) { if (v[NFLX+i] < 0.0) v[NFLX+i] = 0.0; if (v[NFLX+i] > 1.0) v[NFLX+i] = 1.0; CoolCoeffs.dmuN_dX[i] = elem_mass[elem_part[i]]*elem_ab[elem_part[i]]; CoolCoeffs.dmuD_dX[i] = elem_ab[elem_part[i]] *rad_rec_z[i]; mmw1 += CoolCoeffs.dmuN_dX[i]*v[NFLX+i]; /* Numerator part of mu */ mmw2 += CoolCoeffs.dmuD_dX[i]*v[NFLX+i]; /* Denominator part of mu */ } /* -- Add contributions from ionized H -- */ CoolCoeffs.dmuN_dX[0] += -elem_mass[0]*elem_ab[el_H]; CoolCoeffs.dmuD_dX[0] += -2.0*elem_ab[el_H]; mmw1 += elem_mass[0]*elem_ab[el_H]*(1.0 - v[X_HI]); mmw2 += elem_ab[el_H]*(1.0 - v[X_HI])*2.; CoolCoeffs.muN = mmw1; CoolCoeffs.muD = mmw2; if (mmw1 != mmw1) { print(">>> Error! MMW1 NaN! %ld\n",g_stepNumber); for (i = 0; i < NIONS; i++) { print ("%d %10.4e\n",i,v[NFLX+i]); } QUIT_PLUTO(1); } if (mmw2 != mmw2) { print(">>> Error! MMW2 NaN!\n"); for (i = 0; i < NIONS; i++) { print ("%d %10.4e\n",i,v[NFLX+i]); } QUIT_PLUTO(1); } mu = mmw1/mmw2; #endif return mu; }