void hll_flux(double primL[], double primR[], double F[], int nq, double x, double w, struct parList *pars) { double sL, sR, sC; double U[nq], UL[nq], UR[nq], FL[nq], FR[nq]; prim2cons(primL, UL, x, 1.0, pars); prim2cons(primR, UR, x, 1.0, pars); flux(primL, FL, x, pars); flux(primR, FR, x, pars); wave_speeds(primL, primR, &sL, &sR, &sC, x, pars); int q; if(sL > w) for(q=0; q<nq; q++) { U[q] = UL[q]; F[q] = FL[q]; } else if(sR < w) for(q=0; q<nq; q++) { U[q] = UR[q]; F[q] = FR[q]; } else for(q=0; q<nq; q++) { U[q] = (sR*UR[q] - sL*UL[q] - (FR[q] - FL[q])) / (sR - sL); F[q] = (sR*FL[q] - sL*FR[q] + sL*sR*(UR[q] - UL[q])) / (sR - sL); } for(q=0; q<nq; q++) F[q] -= w*U[q]; }
void hllc_flux(double primL[], double primR[], double F[], int nq, double x, double w, struct parList *pars) { double sL, sR, sC; double U[nq], UL[nq], UR[nq], FL[nq], FR[nq]; prim2cons(primL, UL, x, 1.0, pars); prim2cons(primR, UR, x, 1.0, pars); flux(primL, FL, x, pars); flux(primR, FR, x, pars); wave_speeds(primL, primR, &sL, &sR, &sC, x, pars); int q; if(sL > w) for(q=0; q<nq; q++) { U[q] = UL[q]; F[q] = FL[q]; } else if(sR < w) for(q=0; q<nq; q++) { U[q] = UR[q]; F[q] = FR[q]; } else if(sL <= w && w <= sC) { Ustar(primL, U, sL, sC, x, pars); for(q=0; q<nq; q++) F[q] = FL[q] + sL*(U[q]-UL[q]); } else if(sC < w && w <= sR) { Ustar(primR, U, sR, sC, x, pars); for(q=0; q<nq; q++) F[q] = FR[q] + sR*(U[q]-UR[q]); } else { printf("ERROR: wave speeds in HLLC. sL=%.12lg, s=%.12lg, sR=%.12lg\n", sL, sC, sR); } for(q=0; q<nq; q++) F[q] -= w*U[q]; }
void lax_friedrichs_flux(double primL[], double primR[], double F[], int nq, double x, double w, struct parList *pars) { double sL, sR, sC, s; double U[nq], UL[nq], UR[nq], FL[nq], FR[nq]; prim2cons(primL, UL, x, 1.0, pars); prim2cons(primR, UR, x, 1.0, pars); flux(primL, FL, x, pars); flux(primR, FR, x, pars); wave_speeds(primL, primR, &sL, &sR, &sC, x, pars); s = fabs(sL)>fabs(sR) ? fabs(sL) : fabs(sR); int q; for(q=0; q<nq; q++) U[q] = 0.5*(UL[q] + UR[q] - (FR[q] - FL[q])/s); for(q=0; q<nq; q++) F[q] = 0.5*(FL[q] + FR[q] - s*(UR[q] - UL[q])) - w * U[q]; }
void setupCells( struct domain * theDomain ){ int restart_flag = theDomain->theParList.restart_flag; if( restart_flag ) restart( theDomain ); calc_dp( theDomain ); int i,j,k; struct cell ** theCells = theDomain->theCells; int Nr = theDomain->Nr; int Nz = theDomain->Nz; int * Np = theDomain->Np; int Npl = theDomain->Npl; double * r_jph = theDomain->r_jph; double * z_kph = theDomain->z_kph; int atmos = theDomain->theParList.include_atmos; for( j=0 ; j<Nr ; ++j ){ for( k=0 ; k<Nz ; ++k ){ int jk = j+Nr*k; for( i=0 ; i<Np[jk] ; ++i ){ struct cell * c = &(theCells[jk][i]); double phip = c->piph; double phim = phip-c->dphi; c->wiph = 0.0; double xp[3] = {r_jph[j ],phip,z_kph[k ]}; double xm[3] = {r_jph[j-1],phim,z_kph[k-1]}; double r = get_moment_arm( xp , xm ); double dV = get_dV( xp , xm ); double phi = c->piph-.5*c->dphi; double x[3] = { r , phi , .5*(z_kph[k]+z_kph[k-1])}; if( !restart_flag ) initial( c->prim , x ); subtract_omega( c->prim ); if( atmos ){ int p; for( p=0 ; p<Npl ; ++p ){ double gam = theDomain->theParList.Adiabatic_Index; adjust_gas( theDomain->thePlanets+p , x , c->prim , gam ); } } prim2cons( c->prim , c->cons , x , dV ); cons2prim( c->cons , c->prim , x , dV ); } } } set_wcell( theDomain ); }
void boundary_r( struct domain * theDomain ){ int Nt = theDomain->Nt; int Np = theDomain->Np; int * Nr = theDomain->Nr; double * t_jph = theDomain->t_jph; double * p_kph = theDomain->p_kph; struct cell ** theCells = theDomain->theCells; double t = theDomain->t; double vw = 1.0; double Amp = theDomain->theParList.Explosion_Energy; double wavenumber = theDomain->theParList.Gam_0; int j,k; for( j=0 ; j<Nt ; ++j ){ double tp = t_jph[j]; double tm = t_jph[j-1]; for( k=0 ; k<Np ; ++k ){ double pp = p_kph[k]; double pm = p_kph[k-1]; int jk = j+Nt*k; struct cell * c1 = &(theCells[jk][0]); struct cell * c2 = &(theCells[jk][Nr[jk]-1]); struct cell * c2_c = &(theCells[jk][Nr[jk]-2]); double r1 = c1->riph; double x[3] = {r1,.5*(tp+tm),.5*(pp+pm)}; double xp[3] = {r1, tp , pp}; double xm[3] = {0.0,tm,pm}; double dV = get_dV( xp , xm ); initial( c1->prim , x ); double phase1 = .5*(tp+tm); double phase2 = vw*t/r1; double k_p = wavenumber; c1->prim[RHO] *= 1. + Amp*sin(2.5*k_p*phase1)*sin(k_p*phase2); c1->prim[RHO] *= 1. + 0.*.5*2.*((double)rand()/(double)RAND_MAX-.5); prim2cons( c1->prim , c1->cons , r1 , dV ); int q; for( q=0 ; q<NUM_Q ; ++q ){ c1->RKcons[q] = c1->cons[q]; c2->prim[q] = c2_c->prim[q]; } } } }
void calc_cons(struct grid *g, struct parList *pars) { int i, j; int nx1 = g->nx1; int nx2 = g->nx2; int d1 = g->d1; int d2 = g->d2; for(i=0; i<nx1; i++) for(j=0; j<nx2; j++) { double xm[2] = {g->x1[i], g->x2[j]}; double xp[2] = {g->x1[i+1], g->x2[j+1]}; double x[2]; geom_CM(xm,xp,x); double dV = geom_dV(xm,xp); prim2cons(&(g->prim[d1*i+d2*j]), &(g->cons[d1*i+d2*j]), x, dV, pars); } }
void solve_riemann( double * primL , double * primR , double * consL , double * consR , double * gradL , double * gradR , double * x , double * n , double w , double dAdt , int dim , double * E1_riemann , double * B1_riemann , double * E2_riemann , double * B2_riemann ){ int q; double r = x[0]; double Flux[NUM_Q]; double Ustr[NUM_Q]; for( q=0 ; q<NUM_Q ; ++q ){ Flux[q] = 0.0; Ustr[q] = 0.0; } if( riemann_solver == _HLL_ || riemann_solver == _HLLC_ ){ double Sl,Sr,Ss; double Bpack[5]; vel( primL , primR , &Sl , &Sr , &Ss , n , x , Bpack ); if( w < Sl ){ flux( primL , Flux , x , n ); prim2cons( primL , Ustr , x , 1.0 ); }else if( w > Sr ){ flux( primR , Flux , x , n ); prim2cons( primR , Ustr , x , 1.0 ); }else{ if( riemann_solver == _HLL_ ){ double Fl[NUM_Q]; double Fr[NUM_Q]; double Ul[NUM_Q]; double Ur[NUM_Q]; double aL = Sr; double aR = -Sl; prim2cons( primL , Ul , x , 1.0 ); prim2cons( primR , Ur , x , 1.0 ); flux( primL , Fl , x , n ); flux( primR , Fr , x , n ); for( q=0 ; q<NUM_Q ; ++q ){ Flux[q] = ( aL*Fl[q] + aR*Fr[q] + aL*aR*( Ul[q] - Ur[q] ) )/( aL + aR ); Ustr[q] = ( aR*Ul[q] + aL*Ur[q] + Fl[q] - Fr[q] )/( aL + aR ); } }else{ double Uk[NUM_Q]; double Fk[NUM_Q]; if( w < Ss ){ prim2cons( primL , Uk , x , 1.0 ); getUstar( primL , Ustr , x , Sl , Ss , n , Bpack ); flux( primL , Fk , x , n ); for( q=0 ; q<NUM_Q ; ++q ){ Flux[q] = Fk[q] + Sl*( Ustr[q] - Uk[q] ); } }else{ prim2cons( primR , Uk , x , 1.0 ); getUstar( primR , Ustr , x , Sr , Ss , n , Bpack ); flux( primR , Fk , x , n ); for( q=0 ; q<NUM_Q ; ++q ){ Flux[q] = Fk[q] + Sr*( Ustr[q] - Uk[q] ); } } } } }else{ get_Ustar_HLLD( w , primL , primR , Flux , Ustr , r , n ); } if( visc_flag ){ double vFlux[NUM_Q]; double prim[NUM_Q]; double gprim[NUM_Q]; for( q=0 ; q<NUM_Q ; ++q ){ prim[q] = .5*(primL[q]+primR[q]); gprim[q] = .5*(gradL[q]+gradR[q]); if( dim==0 ) gprim[q] /= r; vFlux[q] = 0.0; } visc_flux( prim , gprim , vFlux , x , n ); for( q=0 ; q<NUM_Q ; ++q ) Flux[q] += vFlux[q]; } for( q=0 ; q<NUM_Q ; ++q ){ consL[q] -= (Flux[q] - w*Ustr[q])*dAdt; consR[q] += (Flux[q] - w*Ustr[q])*dAdt; } flux_to_E( Flux , Ustr , x , E1_riemann , B1_riemann , E2_riemann , B2_riemann , dim ); /* if( use_B_fields && NUM_Q > BZZ ){ if( dim==0 ){ *E1_riemann = Flux[BRR]*r; //Ez *B1_riemann = Ustr[BRR]*r*r; // r*Br *E2_riemann = Flux[BZZ]; //Er *B2_riemann = Ustr[BZZ]*r; //-r*Bz }else if( dim==1 ){ *E1_riemann = -Flux[BPP]*r; //Ez *B1_riemann = Ustr[BRR]*r*r; // r*Br *E2_riemann = 1.0*Flux[BZZ]; //Ephi }else{ *E1_riemann = -Flux[BPP]*r; //Er *B1_riemann = Ustr[BZZ]*r; //-r*Bz *E2_riemann = 1.0*-Flux[BRR]*r; //Ephi } } */ }