Scrambler  1
RadShock/problem.f90
Go to the documentation of this file.
00001 !#########################################################################
00002 !               
00003 !    Copyright (C) 2003-2012 Department of Physics and Astronomy,
00004 !                            University of Rochester,
00005 !                            Rochester, NY
00006 !
00007 !    problem.f90 of module RadShock is part of AstroBEAR.
00008 !
00009 !    AstroBEAR is free software: you can redistribute it and/or modify    
00010 !    it under the terms of the GNU General Public License as published by 
00011 !    the Free Software Foundation, either version 3 of the License, or    
00012 !    (at your option) any later version.
00013 !
00014 !    AstroBEAR is distributed in the hope that it will be useful, 
00015 !    but WITHOUT ANY WARRANTY; without even the implied warranty of
00016 !    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00017 !    GNU General Public License for more details.
00018 !
00019 !    You should have received a copy of the GNU General Public License
00020 !    along with AstroBEAR.  If not, see <http://www.gnu.org/licenses/>.
00021 !
00022 !#########################################################################
00025 
00028 
00032 
00035 MODULE Problem
00036 
00037   USE DataDeclarations
00038   USE PhysicsDeclarations
00039   USE GlobalDeclarations
00040   USE SourceDeclarations
00041   USE CoolingSrc
00042   USE NEQCoolingSrc
00043   USE Emissions
00044   USE Projections
00045   USE Fields
00046   USE Cameras
00047   USE ProcessingDeclarations
00048   IMPLICIT NONE
00049   SAVE
00050   
00051   PUBLIC ProblemModuleInit, ProblemGridInit, &
00052      ProblemBeforeStep, ProblemAfterStep, ProblemSetErrFlag, ProblemBeforeGlobalStep
00053   
00054   TYPE(CoolingDef),POINTER :: coolingobj  
00055   INTEGER :: iCooling
00056   REAL(KIND=qPREC) :: n_amb, v_amb, T_amb, By_amb, Xh_amb, Hefrac, alpha, beta
00057   REAL(KIND=qPREC) :: rho_amb, p_amb, rho_ps, v_ps, p_ps, By_ps, mu
00058 
00059   NAMELIST /ProblemData/ n_amb, v_amb, T_amb, By_amb, Xh_amb, Hefrac, alpha, beta, iCooling
00060 
00061 CONTAINS
00062 
00063   SUBROUTINE ProblemModuleInit()
00064     INTEGER :: j
00065     REAL(KIND=qPREC) :: cs, M, a, b, c, d
00066     COMPLEX(8) :: MyRoots(3)
00067     REAL(KIND=qPREC) :: ValidRoots(3)
00068     TYPE(ProjectionDef), POINTER :: Projection
00069 
00070     OPEN(UNIT=Problem_Data_Handle, FILE='problem.data', STATUS="OLD")
00071     READ(Problem_Data_Handle, NML=ProblemData)
00072     CLOSE(Problem_Data_Handle)
00073 
00074     IF (.NOT. lRestart) THEN
00075        CALL CreateCoolingObject(coolingobj)
00076     ELSE
00077        coolingobj => firstcoolingobj
00078     END IF
00079 
00080     coolingobj%iCooling = iCooling; coolingobj%floortemp = 1d3; coolingobj%mintemp = 1d-2;
00081 
00082     coolingobj%alpha = alpha; coolingobj%beta = beta
00083 
00084     lHe = .TRUE.; !lHeII = .TRUE.
00085 
00086     ! solve for initial postshock conditions ("ps" values)
00087     ! v_ps can be found using ambient mach number M and the shock jump equation in hydro case
00088     ! for MHD, v_ps is found via a more complicated cubic equation
00089     ! then rho_ps and p_ps can be found from mass flux and momentum flux conservation
00090     rho_amb = hMass*n_amb; p_amb = n_amb*Boltzmann*T_amb
00091 
00092     IF(lMHD) THEN
00093        ValidRoots = 2d0*v_amb
00094        a = rho_amb*v_amb*(half-gamma15)
00095        b = gamma15*(rho_amb*v_amb**2d0 + p_amb + By_amb**2d0/8d0/Pi)
00096        c = -(half*rho_amb*v_amb**3d0 + gamma15*p_amb*v_amb + By_amb**2d0*v_amb/4d0/Pi)
00097        d = By_amb**2d0*v_amb**2d0/4d0/Pi*(1d0-gamma15/2d0)
00098        MyRoots = CubicRoots(a,b,c,d)
00099        DO j=1, 3
00100           IF(AIMAG(MyRoots(j)) == 0d0 .AND. REAL(MyRoots(j)) > 0d0) THEN
00101              ValidRoots(j) = REAL(MyRoots(j))
00102           END IF
00103        END DO
00104        v_ps = MINVAL(ValidRoots)
00105        rho_ps = rho_amb*v_amb/v_ps
00106        By_ps = By_amb*v_amb/v_ps
00107        p_ps = rho_amb*v_amb**2d0 + p_amb + By_amb**2d0/(8d0*Pi) - rho_ps*v_ps**2d0 - By_ps**2d0/(8d0*Pi)
00108     ELSE
00109        cs = SQRT(gamma*p_amb/rho_amb)
00110        M = v_amb/cs
00111        v_ps = v_amb*(gamma1*M**2d0 + 2d0)/((gamma+1d0)*M**2d0)
00112        rho_ps = rho_amb*v_amb/v_ps
00113        p_ps = rho_amb*v_amb**2d0 + p_amb  - rho_ps*v_ps**2d0
00114     END IF
00115 
00116     CALL CreateProjection(Projection)
00117     Projection%Field%id = SII_6716_Field
00118     Projection%Field%component = GASCOMP
00119     Projection%Field%name = 'SII_6716'
00120     Projection%dim = 3d0
00121 
00122 !    CALL AddDiagnosticVar(MPI_ID_FIELD, 'MPI_ID')
00123 !    CALL AddDiagnosticVar(ChildMask_Field, 'ChildMask')
00124 !    CALL AddDiagnosticVar(ErrFlag_Field, 'ErrFlags')
00125    
00126 !    CALL AddDiagnosticVar(SII_6716_FIELD, 'SII_6716')
00127 !     CALL StoreEmiss(iSII_6731)
00128 !     CALL StoreEmiss(iOI)         ! stores 6300 line
00129 !     CALL StoreEmiss(iNII)        ! stores 6583 line
00130 !     CALL StoreEmiss(iHalpha)
00131 
00132   END SUBROUTINE ProblemModuleInit
00133 
00134   SUBROUTINE ProblemGridInit(Info)
00135 
00136     TYPE(InfoDef) :: Info
00137     REAL(KIND=qPREC), POINTER, DIMENSION(:,:,:,:) :: q
00138     INTEGER :: rmbc, mx, i
00139     REAL(KIND=qPREC) :: x, Lx, xlower, shock_pos, dx, xold
00140     REAL(KIND=qPREC) :: p_cool, v_cool, rho_cool, By_cool, Eflux_cool, HIIflux_cool, HeIIflux_cool, Heflux_cool
00141 
00142     q => Info%q
00143     q = 0d0
00144 
00145     rmbc=levels(Info%level)%gmbc(levels(Info%level)%step)
00146     mx=Info%mX(1); dx=levels(Info%level)%dX*lscale
00147     xlower=Info%xbounds(1,1)*lscale; Lx=Gxbounds(1,2)-Gxbounds(1,1)
00148     shock_pos = 1d-1*Lx*lscale
00149     xold = shock_pos
00150     Eflux_cool = half*rho_ps*v_ps**3d0 + gamma15*p_ps*v_ps
00151     IF(lMHD) Eflux_cool = Eflux_cool + By_ps**2d0*v_ps/4d0/Pi
00152     HIIflux_cool = Xh_amb*rho_ps*v_ps
00153     HeIIflux_cool = 0d0
00154     Heflux_cool = Hefrac*(muHe/muH)*rho_ps*v_ps
00155 
00156     DO i=1-rmbc, mx+rmbc
00157        x=(xlower + (REAL(i) - half) * dx)
00158 
00159        IF(x < shock_pos) THEN
00160           ! ambient preshock region ("amb" values)
00161           q(i,:,:,1) = n_amb / nScale
00162           q(i,:,:,2) = q(i,:,:,1) * v_amb / VelScale
00163           q(i,:,:,iE) = gamma7*p_amb/pScale + half*q(i,:,:,1)*(v_amb/VelScale)**2d0
00164           IF(lMHD) THEN
00165              q(i,:,:,iBy) = By_amb/BScale
00166              q(i,:,:,iE) = q(i,:,:,iE) + half*q(i,:,:,iBy)**2d0
00167           END IF
00168           
00169           IF(iCooling == NEQCool .OR. iCooling == ZCool) THEN
00170              q(i,:,:,iHe) = Hefrac*q(i,:,:,1)
00171 !             q(i,:,:,iHeII) = 0d0
00172              q(i,:,:,iHII) = Xh_amb*(q(i,:,:,1)-q(i,:,:,iHe))
00173              q(i,:,:,iH) = q(i,:,:,1)-q(i,:,:,iHe)-q(i,:,:,iHII)
00174           END IF
00175 
00176        ELSE
00177           ! postshock region ("cool" values)
00178           ! 4th order Runge-Kutta method is used to find energy flux Eflux_cool
00179           CALL RK4(x-xold, 1d-1*dx, Eflux_cool, HIIflux_cool, HeIIflux_cool, Heflux_cool)
00180           xold = x
00181 
00182           CALL EfluxVars(Eflux_cool, v_cool, rho_cool, By_cool, p_cool)
00183 
00184           q(i,:,:,1) = rho_cool / rScale
00185           q(i,:,:,2) = q(i,:,:,1) * v_cool / VelScale
00186           q(i,:,:,iE) = gamma7*p_cool/pScale + half*q(i,:,:,1)*(v_cool/VelScale)**2d0
00187           IF(lMHD) THEN
00188              q(i,:,:,iBy) = By_cool/BScale
00189              q(i,:,:,iE) = q(i,:,:,iE) + half*q(i,:,:,iBy)**2d0
00190           END IF
00191          
00192           IF(iCooling == NEQCool .OR. iCooling == ZCool) THEN
00193              q(i,:,:,iHII) = HIIflux_cool/v_cool/rScale/muH
00194 !             q(i,:,:,iHeII) = HeIIflux_cool/v_cool/rScale/muHe
00195              q(i,:,:,iHe) = Heflux_cool/v_cool/rScale/muHe
00196              q(i,:,:,iH) = q(i,:,:,1) - q(i,:,:,iHe) - q(i,:,:,iHII) !- q(i,:,:,iHeII)
00197           END IF
00198           
00199        END IF
00200     END DO
00201   END SUBROUTINE ProblemGridInit
00202 
00203     ! This subroutine returns a numerical value for p_cool. It follows the 4th order
00204     ! Runge-Kutta integration method to solve for the pressure profile in the energy flux ODE.
00205     ! RK4(upper bound, pressure at lower bound, desired step size)
00206     SUBROUTINE RK4(dx, h, Eflux, HIIflux, HeIIflux, Heflux)
00207        INTEGER :: i, steps
00208        REAL(KIND=qPREC) :: k1, k2, k3, k4, dx, dh, h, Eflux, ddxEflux
00209        REAL(KIND=qPREC) :: l1, l2, l3, l4, HIIflux, ddxHIIflux
00210        REAL(KIND=qPREC) :: j1, j2, j3, j4, HeIIflux, ddxHeIIflux, Heflux
00211 
00212        ! In case dx/h is not an integer, CEILING function is used, and dh is
00213        ! used to ensure that all of dx is divided into evenly distributed steps.
00214        steps = CEILING(dx/h)
00215        dh = dx/REAL(steps)
00216 
00217        DO i=1, steps
00218           CALL Derivs(Eflux, ddxEflux, HIIflux, ddxHIIflux, HeIIflux, ddxHeIIflux, Heflux)
00219           k1 = dh*ddxEflux ; l1 = dh*ddxHIIflux; j1 = dh*ddxHeIIflux
00220           CALL Derivs(Eflux + half*k1, ddxEflux, HIIflux + half*l1, ddxHIIflux, HeIIflux + half*j1, ddxHeIIflux, Heflux - half*j1)
00221           k2 = dh*ddxEflux ; l2 = dh*ddxHIIflux; j2 = dh*ddxHeIIflux
00222           CALL Derivs(Eflux + half*k2, ddxEflux, HIIflux + half*l2, ddxHIIflux, HeIIflux + half*j2, ddxHeIIflux, Heflux - half*j2)
00223           k3 = dh*ddxEflux ; l3 = dh*ddxHIIflux; j3 = dh*ddxHeIIflux
00224           CALL Derivs(Eflux + k3, ddxEflux, HIIflux + l3, ddxHIIflux, HeIIflux + j3, ddxHeIIflux, Heflux - j3)
00225           k4 = dh*ddxEflux ; l4 = dh*ddxHIIflux; j4 = dh*ddxHeIIflux
00226        
00227           Eflux = Eflux + k1/6d0 + k2/3d0 + k3/3d0 + k4/6d0
00228           HIIflux = HIIflux + l1/6d0 + l2/3d0 + l3/3d0 + l4/6d0
00229           HeIIflux = HeIIflux + j1/6d0 + j2/3d0 + j3/3d0 + j4/6d0
00230           Heflux = Heflux - HeIIflux
00231        END DO
00232        
00233     END SUBROUTINE RK4
00234 
00235     ! This function provides values for the right hand side of the equation in the Runge-Kutta
00236     ! method in order to find the change in energy flux.
00237     SUBROUTINE Derivs(Eflux, ddxEflux, HIIflux, ddxHIIflux, HeIIflux, ddxHeIIflux, Heflux)
00238        REAL(KIND=qPREC) :: ddxEflux, ddxHIIflux, ddxHeIIflux
00239        REAL(KIND=qPREC) :: p, v, rho, By, T, T_floor, lambda, Eflux, HIIflux, HeIIflux, Heflux, ne, Xh, nH
00240        REAL(KIND=qPREC),DIMENSION(0:nSpeciesHi) :: nvec
00241        REAL(KIND=qPREC), DIMENSION(NrHydroVars) :: qin, f
00242 
00243        CALL EfluxVars(Eflux, v, rho, By, p)
00244 
00245        T = p/Boltzmann/(rho/hMass)     
00246        T_floor = 1d3
00247 
00248        SELECT CASE(iCooling)
00249           CASE(NoCool)
00250              lambda = 0d0
00251           CASE(AnalyticCool)
00252              lambda = (rho/hMass)**2d0*alpha*T**beta
00253              IF(T<=T_floor) lambda = 0d0 
00254          CASE(DMCool)
00255              lambda = (rho/hMass)**2d0*DMCoolingRate(T)
00256              IF(T<=T_floor) lambda = 0d0 
00257           CASE(NEQCool)
00258              f = 0d0
00259              qin = 0d0
00260              qin(1) = rho/rScale
00261              qin(iHII) = HIIflux/v/rScale/muH
00262 !             qin(iHeII) = HeIIflux/v/rScale/muHe
00263              qin(iHe) = Heflux/v/rScale/muHe
00264              qin(iH) = qin(1)-qin(iHII)-qin(iHe) !-qin(iHeII)
00265              CALL GetNEQvars(qin, mu, nvec)
00266              T = mu*T
00267              CALL Cool_Derivatives(qin,f,T,nvec)
00268              lambda = f(iE)*pScale/timescale
00269              ddxHIIflux = f(iHII)*rScale/timescale
00270 !             ddxHeIIflux = f(iHeII)*rScale/timescale
00271              ddxHeIIflux = 0d0
00272              IF(T<=T_floor) lambda = 0d0 
00273           CASE(ZCool)
00274              f = 0d0
00275              qin = 0d0
00276              qin(1) = rho/rScale
00277              qin(iHII) = HIIflux/v/rScale/muH
00278 !             qin(iHeII) = HeIIflux/v/rScale/muHe
00279              qin(iHe) = Heflux/v/rScale/muHe
00280              qin(iH) = qin(1) - qin(iHII) - qin(iHe) !- qin(iHeII)
00281              CALL GetNEQvars(qin, mu, nvec)
00282              T = mu*T
00283              CALL Cool_Derivatives(qin,f,T,nvec,ZCool)
00284              CALL GetZvars(nvec,ne,Xh,nH)
00285              IF(ZCoolingRate(ne,T,Xh) == 0d0) THEN
00286                 lambda = f(iE)*pScale/timescale
00287              ELSE
00288                 lambda = (f(iE) + (1d0-Zweight(T))*metal_loss)*pScale/timescale + Zweight(T) * nH**2d0 * ZCoolingRate(ne,T,Xh)
00289              END IF
00290              ddxHIIflux = f(iHII)*rScale*muH/timescale
00291 !             ddxHeIIflux = f(iHeII)*rScale*muHe/timescale
00292              ddxHeIIflux = 0d0
00293              IF(T<=T_floor) lambda = 0d0
00294        END SELECT
00295 
00296        ddxEflux = -lambda
00297 
00298     END SUBROUTINE Derivs
00299 
00300     ! Given the energy flux, returns hydro variables
00301     SUBROUTINE EfluxVars(Eflux, v, rho, By, p)
00302        INTEGER :: j
00303        REAL(KIND=qPREC) :: Eflux, v, rho, By, p, a, b, c, d
00304        COMPLEX(8) :: MyRoots(3)
00305        REAL(KIND=qPREC) :: ValidRoots(3)
00306 
00307        IF(lMHD) THEN
00308           ValidRoots = 2d0*v_amb
00309           a = rho_ps*v_ps*(half-gamma15)
00310           b = gamma15*(rho_ps*v_ps**2d0 + p_ps + By_ps**2d0/8d0/Pi)
00311           c = -Eflux
00312           d = By_ps**2d0*v_ps**2d0/4d0/Pi*(1d0-gamma15/2d0)
00313           MyRoots = CubicRoots(a,b,c,d)
00314           DO j=1, 3
00315              IF(AIMAG(MyRoots(j)) == 0d0 .AND. REAL(MyRoots(j)) > 0d0) THEN
00316                 ValidRoots(j) = REAL(MyRoots(j))
00317              END IF
00318           END DO
00319           v = MINVAL(ValidRoots)
00320           rho = rho_ps*v_ps/v
00321           By = By_ps*v_ps/v
00322           p = rho_ps*v_ps**2d0 + p_ps + By_ps**2d0/(8d0*Pi) - rho*v**2d0 - By**2d0/(8d0*Pi)
00323        ELSE
00324           a = rho_ps*v_ps*(half-gamma15)
00325           b = gamma15*(rho_ps*v_ps**2d0 + p_ps)
00326           c = -Eflux
00327           v = (-b + SQRT(b**2d0 - 4d0*a*c))/(2d0*a)
00328           rho = rho_ps*v_ps/v
00329           By = 0d0
00330           p = rho_ps*v_ps**2d0 + p_ps - rho*v**2d0
00331        END IF
00332     END SUBROUTINE EfluxVars
00333 
00334     ! Gives the roots to the cubic equation ax^3 + bx^2 + cx + d = 0
00335     FUNCTION CubicRoots(a, b, c, d)
00336        COMPLEX(8) :: CubicRoots(3)
00337        REAL(KIND=qPREC) :: a, b, c, d, p, q, DD, y1, y2, y3
00338        REAL(KIND=qPREC) :: phi, temp1, temp2, y2r, y2i, u, v
00339 
00340        ! Step 1: Calculate p and q --------------------------------------------
00341        p  = c/a - b*b/a/a/3d0
00342        q  = (2d0*b*b*b/a/a/a - 9d0*b*c/a/a + 27d0*d/a) / 27d0
00343 
00344        ! Step 2: Calculate DD (discriminant) ----------------------------------
00345        DD = p*p*p/27d0 + q*q/4d0
00346 
00347        ! Step 3: Branch to different algorithms based on DD ------------------
00348 
00349        IF(DD .LT. 0d0) THEN
00350           !       Step 3b:
00351           !       3 real unequal roots -- use the trigonometric formulation
00352           phi = ACOS(-q/2d0/SQRT(ABS(p*p*p)/27d0))
00353           temp1=2d0*SQRT(ABS(p)/3d0)
00354           y1 =  temp1*COS(phi/3d0)
00355           y2 = -temp1*COS((phi+pi)/3d0)
00356           y3 = -temp1*COS((phi-pi)/3d0)
00357           temp1 = b/a/3d0
00358           y2 = y2-temp1
00359           y3 = y3-temp1
00360        ELSE
00361           !       Step 3a:
00362           !       1 real root & 2 conjugate complex roots OR 3 real roots (some are equal)
00363           temp1 = -q/2d0 + SQRT(DD)
00364           temp2 = -q/2d0 - SQRT(DD)
00365           u = ABS(temp1)**(1d0/3d0)
00366           v = ABS(temp2)**(1d0/3d0)
00367           IF(temp1 .LT. 0d0) u = -u
00368           IF(temp2 .LT. 0d0) v = -v
00369           y1  = u + v
00370           y2r = -(u+v)/2d0
00371           y2i =  (u-v)*SQRT(3d0)/2d0
00372           temp1 = b/a/3d0
00373           y2r = y2r-temp1
00374        END IF
00375 
00376         ! Step 4: Final transformation -----------------------------------------
00377 
00378         y1 = y1-temp1
00379 
00380         ! Assign answers -------------------------------------------------------
00381         IF(DD .LT. 0d0) THEN
00382            CubicRoots(1) = CMPLX( y1,  0d0)
00383            CubicRoots(2) = CMPLX( y2,  0d0)
00384            CubicRoots(3) = CMPLX( y3,  0d0)
00385         ELSE IF(DD .EQ. 0d0) THEN
00386            CubicRoots(1) = CMPLX( y1,  0d0)
00387            CubicRoots(2) = CMPLX(y2r,  0d0)
00388            CubicRoots(3) = CMPLX(y2r,  0d0)
00389         ELSE
00390            CubicRoots(1) = CMPLX( y1,  0d0)
00391            CubicRoots(2) = CMPLX(y2r, y2i)
00392            CubicRoots(3) = CMPLX(y2r,-y2i)
00393         END IF
00394 
00395     END FUNCTION CubicRoots
00396 
00397   ! Place any pre-processing operations here
00398   SUBROUTINE ProblemBeforeStep(Info)
00399     TYPE(InfoDef) :: Info
00400   END SUBROUTINE ProblemBeforeStep
00401 
00402   ! Place any post-processing operations here
00403   SUBROUTINE ProblemAfterStep(Info)
00404     TYPE(InfoDef) :: Info
00405   END SUBROUTINE ProblemAfterStep
00406 
00407   ! Can be used to set additional refinement
00408   SUBROUTINE ProblemSetErrFlag(Info)
00409     TYPE(InfoDef) :: Info
00410   END SUBROUTINE ProblemSetErrFlag
00411 
00412   SUBROUTINE ProblemBeforeGlobalStep(n)
00413      INTEGER :: n
00414   END SUBROUTINE ProblemBeforeGlobalStep
00415 
00416 END MODULE Problem
 All Classes Files Functions Variables