c subroutine ifelse(ifdouble,1,d_,)do_ppmde0_1dc ^ ( xl, rho, p, u, e_total, ceul, rhonu, ^ e_totnu, unu, ^ dt, smlrho, smallu, smalle, courmx, ^ n_real_zones_in, nbdy_set) c c ************************************************************************** c * THIS CODE IS DO_PPMDE0_1D. * c * THIS CODE WAS WRITTEN BY: WOODWARD RESEARCH GROUP. (c) 1998 * c * ALL RIGHTS RESERVED. VOID WERE PROHIBITED * c ************************************************************************** c c C WHAT This routine performs a 1D direct Eulerian hydrodynamic c update, with artificial diffusion. c c ASSUMPTIONS A uniform, Cartesian grid is assumed. c A general metal) equation of state is used: C p(i) = p00(i) + gamma1(i) * rho(i) * ei(i) c c c p = p00 + gamma1 * rho * ei c This implies that: c ce**2 = gamma1 * (ei + p / rho) c hence: c ce**2 = ((1 + gamma1) * p - p00 ) / rho c Note that the input data should be consistant with the c relation: c ce**2 = dpdrho + ddpdei c c and of course dpdrho must be positive and ddpdei c nonnegative. However, it is not necessary that the c pressure be positive. c C C BOUNDARIES The necessary number of boundary zones to set (NBDY_SET) c is 5, three for interpolation, one for integration, and c Our Fiend Diffuse takes one. It will be further assumed c that at the outset the 2*NBDY_SET boundary zones [often c referred to as fake zones] contain valid values c appropriate to the user's chosen boundary conditions. C This is the only means of signalling boundary conditions C to this routine The user signals his understanding and c acknowledges this by passing to Do_PPM_DE0_1D c the number of boundary zones set on each side. This c routine verifies that sufficient fake zones have been c provided before continuing. If not, an advisory is c displayed and the program terminated. C C DIMENSIONS c input It will be assumed that the dimensions of all zone-averaged c type quantities { e.g. RHO, P, U, RHONU, PNU, etc, ] are C [1-NBDY_SET:N_REAL_ZONES_IN+NBDY_SET]. C The dimensions of all [and I mean ALL] interface type C quantities [e.g. xl, dvoll, dmassl, pavl, uavl, anything C referenced to and interface, usually the lefthand zone c interface] are C [1-NBDY_SET:N_REAL_ZONES_IN+NBDY_SET+1]. C Note the "+1" please. There may be a quiz later. C C output At the completion of this entire routine, new zone-averaged C quantities [e.g. rhonu, pnu, unu] will have updated values C on the range C [1:N_REAL_ZONES_IN], C and any returned interface type thingies will have valid C values on the range [There are none at this time]: C [1:N_REAL_ZONES_IN+1]. c C c C DATA C input Input arrays are: C C xl: the locations of the lefthand zone interfaces, c defining the computational grid. Please realize c that a numerical grid of n zones would require n+1 xl c values (interface locations.) C rho: The zone-averaged densities C p: The zone-averaged pressures C u: The zone-averaged velocities C e_total: the zone-averaged total energy C ceul: A zone-averaged eulerian sound speed c c A general metal) equation of state is used: C p(i) = p00(i) + gamma1(i) * rho(i) * ei(i) c p = p00 + gamma1 * rho * ei c This implies that: c ce**2 = gamma1 * (ei + p / rho) c hence: c ce**2 = ((1 + gamma1) * p - p00 ) / rho c Note that the input data should be consistant with the c relation: c ce**2 = dpdrho + ddpdei c c and of course dpdrho must be positive and ddpdei c nonnegative. However, it is not necessary that the c pressure be positive. c c c C Input scalers are: C C dt: the timestep. C smlrho: a trivial value for density C smallu: a trivial value for velocity C smalle: a trivial value for nergy C the above 4 values are generally about C 1.0e-06 for single precision and C 1.0d-08 for dourble precision. C n_real_zones_in: The number zones (NOT counting fake zones) to be c processed. Upon completion of the timestep, c this many zones will have valid entries in the c output arrays. C nbdy_set: Number of boundary zones on either end set. This c is a check to ensure sufficient values have been c set by the calling routine. The number of c boundary zones required on either side must be c at least equat to the required values: 5 C C output Output arrays are: C rhonu: The new zone-averaged densities at the new c timestep. C e_totnu: The new zone-averaged total energies at the new c timestep. C unu: The new zone-averaged velocities at the new c timestep. C C Output scaler is: C C courmx: The maximum Courant number found by this call C to this subroutine. C C C C VORSICHT HABEN! C C IT IS ASSUMED THAT APPROPRIATE VALUES OF THE DEPENDENT QUANTITIES C (xl, rho, p ,u) HAVE BEEN PLACED INTO THESE BOUNDARY ZONES C BY THE CALLING ROUTINE TO IMPLIMENT WHATEVER BOUNDARY CONDITIONS C THE USER DESIRED. C C VORSICHT HABEN! C C ALL, REPEAT _*ALL*_ (n_real_zones_in + 2*nbdy_set) C ZONES OF RHO, P, U, and (n_real_zones_in + 2**nbdy_set+1) C ZONES OF XL (note the +1) WILL BE USED AND MUST CONTAIN C VALID VALUES. C C VORSICHT HABEN! C C The user signals his understanding and acknowledges this by C passing to Do_PPM_DE0_1Dthe number of boundary zones C set on each side. Ve haf vays of checking dis. C C C C************************************************************************* C************************************************************************* C************************************************************************* C************************************************************************* C************************************************************************* C************************************************************************* C************************************************************************* C************************************************************************* C C these permit switching btw double and single precision C define(zero,ifelse(ifdouble,1,0.d+00,0.0e+00)) define(one,ifelse(ifdouble,1,1.d+00,1.0e+00)) define(two,ifelse(ifdouble,1,2.d+00,2.0e+00)) define(three,ifelse(ifdouble,1,3.d+00,3.0e+00)) define(four,ifelse(ifdouble,1,4.d+00,4.0e+00)) define(five,ifelse(ifdouble,1,5.d+00,5.0e+00)) define(six,ifelse(ifdouble,1,6.d+00,6.0e+00)) define(seven,ifelse(ifdouble,1,7.d+00,7.0e+00)) define(nine,ifelse(ifdouble,1,9.d+00,9.0e+00)) define(ten,ifelse(ifdouble,1,10.d+00,1.0e+01)) define(twelve,ifelse(ifdouble,1,12.d+00,1.2e+01)) define(ffteen,ifelse(ifdouble,1,15.d+00,1.5e+01)) define(onept1,ifelse(ifdouble,1,1.1d+00,1.1e+00)) define(point9,ifelse(ifdouble,1,.9d+00,9.0e-01)) define(half,ifelse(ifdouble,1,.5d+00,5.0e-01)) define(qtr,ifelse(ifdouble,1,.25d+00,2.5e-01)) define(eighth,ifelse(ifdouble,1,.125d+00,1.25e-01)) define(fifth,ifelse(ifdouble,1,.2d+00,2.0e-01)) define(tenth,ifelse(ifdouble,1,.1d+00,1.0e-01)) define(sxy4th,ifelse(ifdouble,1,.015625d+00,1.5625e-02)) c define(klein,ifelse(ifdouble,1,1.0d-08,1.0e-08)) C************************************************************************* C************************************************************************* C************************************************************************* C************************************************************************* C************************************************************************* C************************************************************************* C************************************************************************* C************************************************************************* C************************************************************************* C ifelse(ifdouble,1,implicit real*8 (a-h,o-z), ) dimension xl(1-nbdy_set:n_real_zones_in + nbdy_set+1), ^ rho(1-nbdy_set:n_real_zones_in + nbdy_set), ^ p(1-nbdy_set:n_real_zones_in + nbdy_set), ^ u(1-nbdy_set:n_real_zones_in + nbdy_set), ^ e_total(1-nbdy_set:n_real_zones_in + nbdy_set), ^ ceul(1-nbdy_set:n_real_zones_in + nbdy_set), ^ rhonu(1-nbdy_set:n_real_zones_in + nbdy_set), ^ e_totnu(1-nbdy_set:n_real_zones_in + nbdy_set), ^ unu(1-nbdy_set:n_real_zones_in + nbdy_set) character*100 msg allocatable :: ^ p00(:), ^ gamma1(:), ^ courno(:), ^ detotl(:), ^ difuse_vel(:), ^ difusl(:), ^ dm(:), ^ dmasll(:), ^ dmasrl(:), ^ dmassl(:), ^ dminv(:), ^ dmnew(:), ^ dmnu(:), ^ dmoml(:), ^ dvol(:), ^ dvolinv(:), ^ dvoll(:), ^ dx(:), ^ dxinv(:), ^ rhoavl(:), ^ stpadb(:), ^ pavl(:), ^ uavl(:), ^ unsmad(:), ^ unsmup(:), ^ niterl(:), ^ cesq(:), ^ clagl(:), ^ rhol(:), ^ pl(:), ^ ul(:), ^ rhor(:), ^ pr(:), ^ ur(:), 2 drho(:), ^ dp(:), ^ du(:), 3 rho6(:), ^ p6(:), ^ u6(:), ^ ullrp(:), ^ pllrp(:), ^ urlrm(:), ^ prlrm(:) small = klein c c We will create and use the following "K" indices as pointers to c the beginning and end of the working array. A count KDO will c also be maintained of the number of zones in the array to c process. Note that all of these shall refer to ZONE counts. c Quantities which depend upon lefthand zone interfaces, will c of course require one additional reference on the right. c the K*1's are appropriate whenever starting from scratch with c interpolation. KSTART1 = 1-NBDY_SET KEND1 = N_REAL_ZONES_IN + NBDY_SET KDO1 = N_REAL_ZONES_IN + 2 * NBDY_SET courmx = 0.0 n_required = 5 if ( nbdy_set .ne. n_required ) then call ppm98_message( ' Do_PPM_DE0_1DC Error. ***',0) call ppm98_message( ^ ' Note. the number of boundary zones set at either end of the ',0) call ppm98_message( ^ ' computational grid do not match the required number.',0) call ppm98_message( ^ ' This routine which performes a single step of a Direct Eulerian',0) call ppm98_message( ^ ' hydrodynamic update (including artificial diffusion) requires',0) write(msg,6666) n_required call ppm98_message( msg,0) 6666 format(1x, i3,' fake zones at each end of the array') call ppm98_message( ^ ' IT IS ASSUMED THAT APPROPRIATE VALUES OF THE DEPENDENT',0) call ppm98_message( ^ ' QUANTITIS (xl, rho, p ,u) HAVE BEEN PLACED INTO THESE',0) call ppm98_message( ^ ' BOUNDARY ZONES BY THE CALLING ROUTINE TO IMPLIMENT ',0) call ppm98_message( ^ ' WHATEVER BOUNDARY CONDITIONS THE USER DESIRED. ALL',0) call ppm98_message( ^ ' (n_real_zones_in + 2*nbdy_set)',0) call ppm98_message( ^ ' ZONES OF RHO, P, U and ',0) call ppm98_message( ^ ' (n_real_zones_in + 2*nbdy_set+1)',0) call ppm98_message( ^ ' ZONES OF XL (note the +1) WILL BE USED AND MUST ',0) call ppm98_message( ^ ' CONTAIN VALID VALUES.',0) write(msg,6667) nbdy_set,n_required call ppm98_message( msg,0) 6667 format(' The calling program passed',i6,' boundary zones but ',i3, ' zones') write(msg,6668) n_real_zones_in call ppm98_message( msg,-1) 6668 format(' are required. There were ',i9,' non-boundary zones passed.') endif allocate ( ^ p00(1-nbdy_set:n_real_zones_in + nbdy_set), ^ gamma1(1-nbdy_set:n_real_zones_in + nbdy_set), ^ courno(1-nbdy_set:n_real_zones_in + nbdy_set), ^ detotl(1-nbdy_set:n_real_zones_in + nbdy_set+1), ^ difuse_vel(1-nbdy_set:n_real_zones_in + nbdy_set), ^ difusl(1-nbdy_set:n_real_zones_in + nbdy_set+1), ^ dm(1-nbdy_set:n_real_zones_in + nbdy_set+1), ^ dmasll(1-nbdy_set:n_real_zones_in + nbdy_set+1), ^ dmasrl(1-nbdy_set:n_real_zones_in + nbdy_set+1), ^ dmassl(1-nbdy_set:n_real_zones_in + nbdy_set+1), ^ dminv(1-nbdy_set:n_real_zones_in + nbdy_set), ^ dmnew(1-nbdy_set:n_real_zones_in + nbdy_set+1), ^ dmnu(1-nbdy_set:n_real_zones_in + nbdy_set+1), ^ dmoml(1-nbdy_set:n_real_zones_in + nbdy_set+1), ^ dvol(1-nbdy_set:n_real_zones_in + nbdy_set), ^ dvolinv(1-nbdy_set:n_real_zones_in + nbdy_set+1), ^ dvoll(1-nbdy_set:n_real_zones_in + nbdy_set), ^ dx(1-nbdy_set:n_real_zones_in + nbdy_set), ^ dxinv(1-nbdy_set:n_real_zones_in + nbdy_set), ^ rhoavl(1-nbdy_set:n_real_zones_in + nbdy_set+1), ^ stpadb(1-nbdy_set:n_real_zones_in + nbdy_set), ^ pavl(1-nbdy_set:n_real_zones_in + nbdy_set+1), ^ uavl(1-nbdy_set:n_real_zones_in + nbdy_set+1), ^ unsmad(1-nbdy_set:n_real_zones_in + nbdy_set), ^ unsmup(1-nbdy_set:n_real_zones_in + nbdy_set), ^ niterl(1-nbdy_set:n_real_zones_in+nbdy_set), ^ cesq(1-nbdy_set:n_real_zones_in+nbdy_set), ^ clagl(1-nbdy_set:n_real_zones_in+nbdy_set), ^ rhol(1-nbdy_set:n_real_zones_in+nbdy_set+1), ^ pl(1-nbdy_set:n_real_zones_in+nbdy_set+1), ^ ul(1-nbdy_set:n_real_zones_in+nbdy_set+1), 1 rhor(1-nbdy_set:n_real_zones_in+nbdy_set+1), ^ pr(1-nbdy_set:n_real_zones_in+nbdy_set+1), ^ ur(1-nbdy_set:n_real_zones_in+nbdy_set+1), 2 drho(1-nbdy_set:n_real_zones_in+nbdy_set), ^ dp(1-nbdy_set:n_real_zones_in+nbdy_set), ^ du(1-nbdy_set:n_real_zones_in+nbdy_set), 3 rho6(1-nbdy_set:n_real_zones_in+nbdy_set), ^ p6(1-nbdy_set:n_real_zones_in+nbdy_set), ^ u6(1-nbdy_set:n_real_zones_in+nbdy_set), ^ ullrp(1-nbdy_set:n_real_zones_in+nbdy_set+1), ^ pllrp(1-nbdy_set:n_real_zones_in+nbdy_set+1), ^ urlrm(1-nbdy_set:n_real_zones_in+nbdy_set+1), 1 prlrm(1-nbdy_set:n_real_zones_in+nbdy_set+1)) N = n_real_zones_in KSTART = 1-nbdy_set KEND = N+nbdy_set KDO = N + 2*nbdy_set c c c compute some essential quantities: c do i=kstart1, kend1 dx(i) = xl(I+1) - xl(i) dxinv(i) = 1.0 / dx(i) dvol(i) = dx(i) dvolinv(i) = dxinv(i) c NOTE these last two are true only in cartesian coords. dm(i) = rho(i) * dvol(i) dminv(i) = 1.0 / dm(i) cesq(i) = ceul(i) * ceul(i) ek = 0.5 * ( u(i)*u(i) ) ei = e_total(i) - ek gamma1(i) = cesq(i) / (ei + p(i)/rho(i) ) p00(i) = p(i) - gamma1(i) * rho(i) * ei enddo c Get the fluxes of the conserved hydrodynamical quantities: c These are mass, longitudinal momentum, and energy in c in the following steps: c c 1) Interpolate, with monotonicity constraints the c pressure and longitudinal velocity by constructing c interpolations of the +- Riemann invarients. c (Returns UNSMUP) c 2) Interpolate, with monotonicity constraints and c discontinuity detection, the density by c constructing interpolations of the 0 Riemann c invarient (streamline, or adiabatic invarient.) c (Returns UNSMAD, STPADB) c 3) Averaging the piecewise parabolic interpolations for c density, velocity, and pressure, find the c averages in the appropriate domains of dependence c for a Riemann problem at each zone interface. c 4) Solve Riemann shock tube problems at each lefthand c zone interface to obtain appropriate time-averaged c values for density, pressure, and longitudinal c momentum at these interfaces. c (Returns RHOAVL, PAVL, UAVL.) c 5) From these time averages, construct interface fluxes c (in the Eulerian sense) of mass, momentum, and c total energy. c (Returns DMASSL, DMOML, and DETOTL.) c 6) Determine the Courant condition (stability) for c each _INTERFACE_. Simply put, in order for c the method to be stable, the various Domains c of dependence for the Riemann problems c cannot exceed the width of either the zone c to the left or the zone to the right. This c make the courant number array COURNO appear c to be an interface-like quantity. c c These individual steps are also available as separate c modules. c nbdy = nbdy_set nzones = n_real_zones_in K = KSTART1 c c c third = one / three forthd = four * third twelth = qtr * third smlfrc = ifelse(ifdouble,1,5.0d-03,5.0e-03) c c c c Now find the spatially-averaged values in the left and right c domain of dependence for the lefthand zone interface for each zone. c c PPM98_EULER0_DOMAIN_AVG calculates the left and right states for c a Riemann problem to each lefthand zone interface in the c following steps: c 1) Interpolate, with monotonicity constraints the c pressure and longitudinal velocity by constructing c interpolations of the +- Riemann invarients. c (Returns UL, UR, DU, U6, PL, PR, DP, P6, UNSMUP) c 2) Interpolate, with monotonicity constraints and c discontinuity detection, the density by c constructing interpolations of the 0 Riemann c invarient (streamline, or adiabatic invarient.) c (Returns UNSMAD, STPADB, RHOL, RHOR, DRHO, RHO6) c 3) Using the piecewise parabolic interpolations for c velocity, and pressure, find the averages in the c appropriate domains of dependence c c c Note: c The average of density over the r0 DOD cannot c be found until the Riemann problem has been c solved. Therefore, the density interplotions c are returned only. nnbdy_set = 4 K = KSTART N = KDO - 2*nnbdy_set call ifelse(ifdouble,1,d_,) ppm98_Euler0_domain_avg ^ ( dx(K), rho(K), p(K), u(K), ceul(K), ^ ullrp(K), pllrp(K), urlrm(K), prlrm(K), clagl(K), ^ rhol(K), rhor(K), drho(K), rho6(K), ul(K), ^ ur(K), du(K), u6(K), pl(K), pr(K), ^ dp(K), p6(K), unsmad(K), stpadb(K), unsmup(K), ^ courno(K), dt, N, nnbdy_set) c c Update the zone counts c kstart2 = kstart+nnbdy_set kend2 = kend-nnbdy_set kdo2 = kdo - 2*nnbdy_set c If one wished to examine the domain averages one could : c c print *,' i, ullrp(i), pllrp(i), urlrm(i), prlrm(i)' c do i=kstart2,kend2+1 c write(6,66) i, ullrp(i), pllrp(i), urlrm(i), prlrm(i) c enddo c print *,' i, rhol(i), rho(i), rhor(i), drho(i), rho6(i)' c do i=kstart2,kend2+1 c write(6,66) i, rhol(i), rho(i), rhor(i), drho(i), rho6(i) c enddo c c c *********************************************************************** c * * c * at this point, modifications to pllrp, prlrm to include body * c * forces, or to ullrp, urlrm to include geometrical effects * c * may be done. (See Colella and Woodward 1984). By modifying * c * the domain averages in this fashon, we may use a single * c * Riemann solver. * c * * c *********************************************************************** c c c Solve the Riemann problems and get the interface fluxes. c c c PPM98_E_RIEMANN uses the domain averages of pressure and c velocity to solve Riemann's shocktube problem at each zone c interface, to provide time-averaged values of pressure c and velocity (returns PAVL and UAVL) at the zone interfaces. c The density pressure interpolations are then used to find c the time-averaged density (returns RHOAVL). From these c values, advected fluxes of volume, mass, momentum, and total c energy are computed. (Returns DVOLL, DMASSL, DMOML, and c DETOT1L) c nnbdy = 1 c KSTART2, KEND2, KDO2 refer to the dod interface values c c we require nnbdy additional zone on the left and right c of the density and pressure interpolations. We will use c these as our points of reference. Hence we adjust these \c c K values KSTART_RIEMANN = KSTART2-nnbdy KEND_RIEMANN = KEND2+nnbdy KDO_RIEMANN = KDO2+2*nnbdy c NITERL is an integer array specifying the number of non-linear iterations c to take at each lefthand zone interface location. In smooth flow, or if c one intends to rely upon the later addition of additional dissapation, no c iterations are necessary. (NITERL(I)=0 for all i). However, situations c may arise inwich the user may wish to take additional iterations. Generally c only in the most extreme situations would more than 1 additional iteration c be required. The results from shock detection may be used to identify c such zones. Here, we rely upon additional dissapation. c K=KSTART_RIEMANN KEND_RIEMANN = KEND2+nnbdy N = KDO_RIEMANN - 2*nnbdy do i=k,k+n+1 niterl(i) = 0 enddo call ifelse(ifdouble,1,d_,) PPM98_E_RIEMANN ^ ( dx(K), ullrp(K), pllrp(K), urlrm(K), prlrm(K), ^ clagl(K), rhol(K), rho(K), rhor(K), drho(K), ^ rho6(K), pl(K), p(K), pr(K), dp(K), ^ p6(K), p00(K), gamma1(K), rhoavl(K), pavl(K), ^ uavl(K), dvoll(K), dmassl(K), dmoml(K), detotl(K), ^ smlrho, dt, niterl(k), n, nnbdy) c c Update the zone counts c kstart3 = kstart2 + nnbdy kend3 = kend2 - nnbdy c c If one wished to examine the fluxes one could: c (Note, these are interface values, hence the +1) c print *,' i, rhoavl(i), pavl(i), uavl(i)' c do i=kstart3,kend3+1 c write(6,66) i, rhoavl(i), pavl(i), uavl(i) c enddo c c print *,' i, dvoll(i), dmassl(i), dmoml(i), detotl(i)' c do i=kstart3,kend3+1 c write(6,66) i, dvoll(i), dmassl(i), dmoml(i), detotl(i) c enddo c c c c c *************************************************************************** * * * c * The conservation laws may now be applied to obtain new zone-averaged * c * values of density, energy, and velocity. At this point, new pressures * c * are not possible because gases from adjacent zones have mixed. The * c * new equation of stae is not known. * c * * c *************************************************************************** c c c update the zone counts. These are the valid counts for the interface c fluxes, dmassl, dvoll, etc. The number of these boundary fluxes which c are for zone THIS program considers fake are (beginning #) - (used #) c This program began with nbdy_set, and nbdy were just used. c nbdy_hydroflux is the number of fake zones with fluxes. nbdy_hydroflux = nbdy_set - nbdy KSTART2 = KSTART1 + nbdy KEND2 = KEND1 - nbdy KDO2 = KDO1 - 2*nbdy c c Thus the fluxes are now valid on the range (kstart2,kend2. To examine c them one could do: c print *,' i, rhoavl(i), pavl(i), uavl(i)' c do i=kstart2,kend2+1 c write(6,666) i, rhoavl(i), pavl(i), uavl(i) c enddo 666 format(1h ,i4,1p,9(1x,e10.3)) 66 format(1h ,i4,1p,9(1x,e10.3)) c the K*2 will be appropriate for working with fluxes (before diffuse) c find the courant number from the above step. When a routine is c called with an array such as COURNO, the old values of COURNO are c overwritten. Thus the Courant number should be determined after each c such routine, or different arrays used. c c To be stable when a flux at xl(i), the left hand edge of zone I, is c computed, the courant numbers in zone I-1 and I need to be <= 1.0 c Thus we examine Courant numbers over the KSTART2-1,KEND2+1 c courmxh = courno(KSTART2 -1) do i=KSTART2 , KEND2+1 courmxh = max (courmxh, courno(i)) enddo c c c PPM98_DEUL_UPDATE_CART upates the zone averages by conserviatvely c differencing the fluxes found above, assuming cartesian c coordinates and no body forces. Non-cartesian geometries c may be treated by including the area terms into the interface c fluxes BEFORE this routine is called. Body forces would c entail and additional term added to the new zone-averages of c velocity and total energy AFTER this routine is called. c c K = KSTART2 nzones = KDO2 call ifelse(ifdouble,1,d_,)ppm98_deul_update_cart ^ ( dx(K), dvol(K), dm(K), rho(K), u(K), ^ e_total(K), dmassl(K), dmoml(K), pavl(K), detotl(K), ^ rhonu(K), unu(K), e_totnu(K), dmnu(K), smlrho, ^ smallp, smallu, smalle, dt, nzones) c no fake zones necessary in the above routine, so the "K" count is c unchanged. The updated zone-averages are valid on the range c kstart2,kend2, and could be examined by: c C print *,' Finished with the hydro step.' C print *,' i, rhonu(i), unu(i), e_totnu(i)' C do i=kstart2,kend2 C write(6,666) i, rhonu(i), unu(i), e_totnu(i) C enddo c C ********************************************************** C * * C * Done with the hydro step, begin diffusion step * C * * C ********************************************************** c zone centerred difusion velocities are required on the range c of zone just updated. Any bondary requirements (offsets) c are built in into the offsets in the calling statement. c c Therefore the correct relative beginning locations is KSTART2 c and the number of zone to process is KDO2 c c BUT REALIZE THAT PPM98_EUL_DIFFUSE_1D OPERATES ON THE INITIAL C DATA. c c c KSTARTDF = KSTART2 KENDDF = KEND2 KDODF = KDO2 c c c Get difusion velocities c c PPM98_DIFUSE0_1D obtains zone centered Eulerian diffusion velocities. c (It's like sausage, you don't really want to know it's done) c K=KSTARTDF nzones = KDODF call ifelse(ifdouble,1,d_,) ppm98_eul_difuse_1d ^ ( dx(k-2), dx(k-1), dx(K), dx(k+1), dx(k+2), ^ u(k-2), u(k-1), u(K), u(k+1), u(k+2), ^ p(k-2), p(K), p(k+2), rho(k-2), rho(K), rho(k+2), ^ ceul(k-2), ceul(K), ceul(k+2), ^ difuse_vel(K), smallu, nzones) c update the zone counts KSTARTDF1 = KSTARTDF KENDDF1 = KENDDF KDODF1 = KDODF c At this point, we should have one extra difuseion velocity c at each end of the working arry. These are ZONE-CENTERED c velocities, from which INTERFACE centerred fluxes will c be constructed. This construction requires this extra fake zone. C Now apply diffusion c c First get the diffusive zone interface fluxes c from the zone centerred difusion velocities nbdy_difuse_fluxes = 1 K = KSTARTDF1 nzones = KDODF1 - 2*nbdy_difuse_fluxes call ifelse(ifdouble,1,d_,) ppm98_difuseflux_hydro_cart ^ ( dx(K), difuse_vel(K), ^ rhonu(K), unu(K), e_totnu(K), ^ dvoll(K), dmasll(K), dmasrl(K), ^ dmassl(K), dmoml(K), detotl(K), ^ courno(K), dt, nzones, ^ nbdy_difuse_fluxes ) KSTARTDF2 = KSTARTDF1 + nbdy_difuse_fluxes KENDDF2 = KENDDF1 - nbdy_difuse_fluxes courdf2 = courno(KSTARTDF2) do i=KSTARTDF2+1,KENDDF2 courdf2 = max (courdf2,courno(i)) enddo c Now, apply the diffusive fluxes c c (Final update:, whew!) c c We occasionally get underflows by allowing absurdly small c velocities to be computed. Therefore, we will set such c velocities to zero below. c ssmalu = small * smallu ssmlu2 = ssmalu * ssmalu c c print *,' i, rhonu(i), unu(i), e_totnu(i), pnu(i)' do 9100 i = KSTARTDF2,KENDDF2 c dmnew(i) = dmnu(i) c this saves a flop or two: if ((dvoll(i) + dvoll(i+1)) .ne. zero) then c dmnew(i) = rhonu(i) * dvol(I) dmnu(i) = dmnew(i) + dmassl(i) - dmassl(i+1) rhonu(i) = dmnu(i) * dvolinv(i) dmnui = one / dmnu(i) unu(i) = (unu(i) * dmnew(i) + dmoml(i) - dmoml(i+1)) * dmnui e_totnu(i) = (e_totnu(i) * dmnew(i) + ^ detotl(i) - detotl(i+1)) * dmnui 9099 continue endif usq = unu(i) * unu(i) c if ((unu(i) * unu(i)) .lt. ssmlu2) unu(i) = 0.0 if ( usq .lt. ssmlu2) then unu(i) = 0.0 usq = 0.0 endif c write(6,666) i, rhonu(i), unu(i), e_totnu(i) 9100 continue c c The final Courant number is the maximum of the c hydrodynamics Courant number and the diffusive c Courant number. c courmx = max (courmxh, courdf2 ) c deallocate ( ^ p00, ^ gamma1, ^ courno, ^ detotl, ^ difuse_vel, ^ difusl, ^ dm, ^ dmasll, ^ dmasrl, ^ dmassl, ^ dminv, ^ dmnew, ^ dmnu, ^ dmoml, ^ dvol, ^ dvolinv, ^ dvoll, ^ dx, ^ dxinv, ^ rhoavl, ^ stpadb, ^ pavl, ^ uavl, ^ unsmad, ^ unsmup, ^ niterl, ^ cesq, ^ clagl, ^ rhol, ^ pl, ^ ul, ^ rhor, ^ pr, ^ ur, 2 drho, ^ dp, ^ du, 3 rho6, ^ p6, ^ u6, ^ ullrp, ^ pllrp, ^ urlrm, ^ prlrm) return end