Scrambler  1
cooling.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 !    cooling.f90 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 
00029 !! @ingroup Source
00030 
00032 MODULE CoolingSrc
00033 
00034   USE DataDeclarations
00035   USE PhysicsDeclarations
00036   USE EOS
00037   USE SourceDeclarations
00038   USE NEQCoolingSrc
00039 
00040   IMPLICIT NONE
00041   PRIVATE
00042   PUBLIC Cooling,CreateCoolingObject,DestroyCoolingObject,InitCoolingTracers, InitCoolingElliptics, CoolingCheck, CoolingIO, CoolingFinalizeInit, Cooling_CountObjects
00043   PUBLIC Cooling_ReadObjectFromChombo, Cooling_InitChomboDatasets, Cooling_WriteObjectToChombo, GetIICoolEqTemp, InitIICool, IICoolingRate, GetCoolingStrength, GetIICoolEqDensity, InitZCool, GetZvars, ZCoolingRate, Zweight
00044   TYPE,PUBLIC :: CoolingDef
00045      ! analytic cooling alpha, beta
00046      REAL(KIND=qPrec) :: alpha,beta
00047      INTEGER :: iCooling    ! see parameters for iCooling below
00048      LOGICAL :: lPrimitive  ! primitive or conservative form of source term
00049      REAL(KIND=qPrec) :: FloorTemp ! do not apply cooling below this temperature
00050      REAL(KIND=qPrec) :: MinTemp   ! do not allow things to cool below this temperature
00051      !
00052      ! parameters to specify region to apply source term, *** spatially driven ***
00053      REAL(KIND=xPrec) :: xlower(3)=0d0,xupper(3)=0d0  ! rectangular region
00054      REAL(KIND=xPrec) :: pos(3),radius=0d0            ! spherical region
00055      ! parameters to specify range to apply source term, *** variable driven ***
00056      INTEGER,DIMENSION(:),ALLOCATABLE :: var
00057      REAL(KIND=qPrec),DIMENSION(:,:),ALLOCATABLE :: varRange
00058      ! scaling to computational units
00059      REAL(KIND=qPrec) :: ScaleCool
00060      ! in locations where two source terms overlap,
00061      ! one should have priority (more than two will not work right)
00062      LOGICAL :: lHasPriority
00063      ! linked-list bookeeping variables
00064      INTEGER :: id
00065      TYPE(CoolingDef),POINTER :: previous,next
00066   END TYPE CoolingDef
00067 
00068   ! Parameters for iCooling
00069   INTEGER,PARAMETER,PUBLIC :: NoCool=0, AnalyticCool=1, DMCool=2, IICool=3, NEQCool=4, ZCool=5
00070 
00071   TYPE(CoolingDef),PUBLIC,POINTER :: firstcoolingobj,lastcoolingobj
00072   INTEGER :: iCoolingObjID=0
00073 
00074   ! Cooling tables
00075   REAL(KIND=qPrec),DIMENSION(:),ALLOCATABLE :: IICoolingTab,lognetab,temptab,logxtab
00076   REAL(KIND=qPrec),DIMENSION(:,:,:),ALLOCATABLE :: ZCoolingTab
00077 
00078   ! Parameters for ZCooling.tab to be read in and used properly
00079   INTEGER :: nDensities, nTemps, nXs
00080   REAL(KIND=qPrec) ::  lognemin, lognemax, tempmin, tempmax, logxmin, logxmax
00081   REAL(KIND=qPREC), PUBLIC :: IICoolPar(1:7)=(/2d-26,1d7,1.184d5,1d3,1.4d-2,9.2d1,0d0/)
00082 CONTAINS
00083 
00084 
00085   ! ==================================================================
00086   ! =                       Cooling IO                               =
00087   ! ==================================================================
00088 
00091   SUBROUTINE CoolingIO(i_IO)
00092     ! Interface declarations
00093     INTEGER,INTENT(IN) :: i_IO
00094     ! Internal declarations
00095     TYPE(CoolingDef),POINTER :: coolingobj
00096 
00097     IF(i_IO==0) THEN
00098        ! output
00099        CALL OutputCooling
00100     ELSE
00101        ! input
00102        CALL InputCooling
00103     END IF
00104   END SUBROUTINE CoolingIO
00105 
00107   SUBROUTINE OutputCooling
00108     ! Internal declarations
00109     TYPE(CoolingDef),POINTER :: coolingobj
00110 
00111     coolingobj=>firstcoolingobj
00112     DO WHILE(ASSOCIATED(coolingobj))
00113        ! Do the output ...
00114 
00115        coolingobj=>coolingobj%next
00116     END DO
00117   END SUBROUTINE OutputCooling
00118 
00120   SUBROUTINE InputCooling
00121     ! Internal declarations
00122     TYPE(CoolingDef),POINTER :: coolingobj
00123 
00124     ! Allocate object, nullify pointers, add to list of cooling objects
00125     CALL CreateCoolingObject(coolingobj)
00126     
00127     ! Do the input ...
00128   END SUBROUTINE InputCooling
00129 
00131   INTEGER FUNCTION Cooling_CountObjects()
00132 
00133       TYPE(CoolingDef), POINTER :: coolingobj
00134       INTEGER :: counter
00135 
00136       coolingobj => firstcoolingobj
00137       counter = 0
00138 
00139       ! Count the number of cooling objects in this list.
00140       DO WHILE (ASSOCIATED(coolingobj))
00141           counter = counter + 1
00142           coolingobj => coolingobj%next
00143       END DO
00144 
00145       Cooling_CountObjects = counter
00146 
00147   END FUNCTION Cooling_CountObjects
00148 
00149   ! ==================================================================
00150   ! =                Main Cooling Section                            =
00151   ! ==================================================================
00152 
00153 
00159   SUBROUTINE Cooling(q,dqdt,x,dx,lform)
00160     ! Interface declarations
00161     REAL(KIND=qPrec) :: q(:)
00162     REAL(KIND=qPrec) :: dqdt(:),x(3),dx
00163     ! Internal declarations
00164     TYPE(CoolingDef),POINTER :: coolingobj
00165     LOGICAL :: lCool
00166     LOGICAL :: lform
00167     coolingobj=>firstcoolingobj
00168     DO WHILE(ASSOCIATED(coolingobj))
00169        ! first determine if cell is within cooling region/regime
00170        !lcool=.FALSE.
00171        !lcool=TestRegion(coolingobj,x,dx)
00172        !lcool=TestRegime(coolingobj,q)
00173 
00174        !IF(lCool) THEN
00175           SELECT CASE(coolingobj%iCooling)
00176           CASE(NoCool)
00177              ! do nothing
00178           CASE(AnalyticCool)
00179              CALL AnalyticCooling(q,dqdt,coolingobj,lform)
00180           CASE(DMCool)
00181              CALL DMCooling(q,dqdt,coolingobj,lform)
00182           CASE(IICool)
00183              CALL IICooling(q,dqdt,coolingobj,lform)
00184           CASE(NEQCool)
00185              CALL NEQCooling(q,dqdt,coolingobj,lform)
00186           CASE(ZCool)
00187              CALL ZCooling(q,dqdt,coolingobj,lform)
00188           CASE DEFAULT
00189           END SELECT
00190        !END IF
00191 
00192        coolingobj=>coolingobj%next
00193     END DO
00194   END SUBROUTINE Cooling
00195 
00196   FUNCTION GetCoolingStrength(q,lform)
00197     REAL(KIND=qPrec) :: q(:)
00198     REAL(KIND=qPrec) :: x(3),dx,GetCoolingStrength
00199     ! Internal declarations
00200     TYPE(CoolingDef),POINTER :: coolingobj
00201     LOGICAL :: lCool
00202     LOGICAL :: lform
00203     coolingobj=>firstcoolingobj
00204     GetCoolingStrength=0d0
00205     DO WHILE(ASSOCIATED(coolingobj))
00206        SELECT CASE(coolingobj%iCooling)
00207        CASE(NoCool)
00208           ! do nothing
00209        CASE(AnalyticCool)
00210           GetCoolingStrength=GetCoolingStrength+AnalyticCoolingStrength(q,coolingobj,lform)
00211        CASE(DMCool)
00212           GetCoolingStrength=GetCoolingStrength+DMCoolingStrength(q,coolingobj,lform)
00213        CASE(IICool)
00214           GetCoolingStrength=GetCoolingStrength+IICoolingStrength(q,coolingobj,lform)
00215        CASE(ZCool)
00216           GetCoolingStrength=GetCoolingStrength+ZCoolingStrength(q,coolingobj,lform)
00217        CASE DEFAULT
00218        END SELECT
00219        coolingobj=>coolingobj%next
00220     END DO
00221 
00222  END FUNCTION GetCoolingStrength
00223 
00224 
00229   LOGICAL FUNCTION TestRegion(coolingobj,x,dx)
00230     TYPE(CoolingDef),POINTER :: coolingobj
00231     REAL(KIND=xPrec) :: x(3),dx
00232     TestRegion=.true.
00233     IF(ALL(coolingobj%xLower==coolingobj%xUpper) .AND. coolingobj%radius==0d0) RETURN
00234   END FUNCTION TestRegion
00235 
00239   LOGICAL FUNCTION TestRegime(coolingobj,q)
00240     TYPE(CoolingDef),POINTER :: coolingobj
00241     REAL(KIND=qPrec)  :: q(:)
00242     TestRegime=.TRUE.
00243     IF(.NOT.(ALLOCATED(coolingobj%var))) RETURN
00244   END FUNCTION TestRegime
00245   
00246 
00247 
00248   ! ==================================================================
00249   ! =               Actual Cooling Source Terms                      =
00250   ! ==================================================================
00251 
00258   SUBROUTINE AnalyticCooling(q,dqdt,coolingobj,lform)
00259     ! Interface declarations
00260     REAL(KIND=qPrec) :: q(:)
00261     REAL(KIND=qPrec),INTENT(INOUT) :: dqdt(:)
00262     TYPE(CoolingDef),POINTER :: coolingobj
00263     LOGICAL :: lform
00264     ! Local declarations
00265 
00266     IF (lform .eqv. PRIMITIVE) THEN      
00267        dqdt(iE) = dqdt(iE)-AnalyticCoolingStrength(q,coolingobj,lform) * gamma1
00268     ELSE
00269        dqdt(iE) = dqdt(iE)-AnalyticCoolingStrength(q,coolingobj,lform)
00270     END IF
00271   END SUBROUTINE AnalyticCooling
00272 
00273   FUNCTION AnalyticCoolingStrength(q,coolingobj,lform)
00274      REAL(KIND=qPREC) :: AnalyticCoolingStrength
00275      REAL(KIND=qPrec) :: q(:)
00276      TYPE(CoolingDef),POINTER :: coolingobj
00277      ! Local declarations
00278      REAL(KIND=qPrec) :: P,T0,dqlocal,alpha,beta,qnew(NrVars),Pnew,Tnew,coolfrac
00279      LOGICAL :: lform
00280 
00281      alpha=coolingobj%alpha
00282      beta=coolingobj%beta
00283      P=MERGE(PrimPress(q),Press(q),lform .eqv. PRIMITIVE)
00284      T0 = P / q(1)
00285      !    write(*,*) alpha, beta, T0, TempScale, coolingobj%floortemp
00286      IF(T0*TempScale<coolingobj%floortemp) THEN
00287         AnalyticCoolingStrength=0d0
00288      ELSE
00289         AnalyticCoolingStrength=q(1)**2 * alpha*T0**beta*coolingobj%ScaleCool
00290      END IF
00291 !     write(*,*) T0, AnalyticCoolingStrength
00292   END FUNCTION AnalyticCoolingStrength
00293 
00294 
00295 
00302   SUBROUTINE DMCooling(q,dqdt,coolingobj,lform)
00303     ! Interface declarations
00304     REAL(KIND=qPrec) :: q(:),dqdt(:)
00305     TYPE(CoolingDef),POINTER :: coolingobj
00306     ! Local declarations
00307     REAL(KIND=qPrec) :: P,T0,dqlocal
00308     LOGICAL :: lform
00309 
00310     IF (lform .eqv. PRIMITIVE) THEN
00311        dqdt(iE) = dqdt(iE) - DMCoolingStrength(q,coolingobj,lform)*gamma1
00312     ELSE
00313        dqdt(iE) = dqdt(iE) - DMCoolingStrength(q,coolingobj,lform)
00314     END IF
00315 
00316   END SUBROUTINE DMCooling
00317 
00318   FUNCTION DMCoolingStrength(q,coolingobj,lform)
00319      REAL(KIND=qPrec) :: q(:),DMCoolingStrength
00320      TYPE(CoolingDef),POINTER :: coolingobj
00321      ! Local declarations
00322      REAL(KIND=qPrec) :: P,T0,dqlocal
00323      LOGICAL :: lform
00324      P=MERGE(PrimPress(q),Press(q),lform .eqv. PRIMITIVE)
00325 
00326      T0 = P/q(1) * TempScale
00327      IF (T0 < max(100d0,coolingobj%floortemp)) THEN
00328         DMCoolingStrength=0d0
00329      ELSE
00330         DMCoolingStrength=q(1)**(2) * DMCoolingRate(T0) * coolingobj%ScaleCool
00331      END IF
00332   END FUNCTION DMCoolingStrength
00333 
00342   SUBROUTINE IICooling(q,dqdt,coolingobj,lform)
00343      ! Interface declarations
00344      REAL(KIND=qPrec) :: q(:),dqdt(:)
00345      TYPE(CoolingDef),POINTER :: coolingobj
00346      LOGICAL :: lform
00347      IF (lform .eqv. PRIMITIVE) THEN
00348         dqdt(iE)=dqdt(iE) + (IIHeatingStrength(q,coolingobj,lform)-IICoolingStrength(q,coolingobj,lform))*gamma1
00349      ELSE
00350         dqdt(iE)=dqdt(iE) + (IIHeatingStrength(q,coolingobj,lform)-IICoolingStrength(q,coolingobj,lform))
00351      END IF
00352   END SUBROUTINE IICooling
00353 
00354   FUNCTION IIHeatingStrength(q,coolingobj, lform)
00355      TYPE(CoolingDef),POINTER :: coolingobj
00356      REAL(KIND=qPREC) :: IIHeatingStrength, q(:)
00357      LOGICAL :: lform
00358      IIHeatingStrength=q(1)*nScale*IICoolPar(1) * coolingobj%ScaleCool
00359   END FUNCTION IIHeatingStrength
00360 
00361   FUNCTION IICoolingStrength(q,coolingobj,lform)
00362      TYPE(CoolingDef),POINTER :: coolingobj
00363      REAL(KIND=qPrec) :: P,T,IICoolingStrength, q(:)
00364      LOGICAL :: lform
00365      P=MERGE(PrimPress(q),Press(q), lform .eqv. PRIMITIVE)
00366      T = P/q(1)*TempScale
00367      IF(T<coolingobj%floortemp) THEN
00368         IICoolingStrength=0d0
00369      ELSE
00370         IICoolingStrength=q(1)**2*nScale**2*IICoolingRate(T)*coolingobj%ScaleCool
00371      END IF
00372   END FUNCTION IICoolingStrength
00373 
00375   FUNCTION GetIICoolEqTemp(rho)
00376 
00377      REAL(KIND=qPrec) :: logtemp,rho, GetIICoolEqTemp
00378      INTEGER :: i
00379      logtemp=4 !(0-8)
00380      !using .9 allows for a range of +- 10 in log space
00381      !If we want 16 digits of accuracy we need -16/log10(.9)=718
00382      DO i=1, 350
00383         IF (rho*IIHeatingRate() - rho**2*IICoolingRate(10d0**logtemp) > 0d0) THEN
00384            logtemp=logtemp+.9d0**i
00385 !           write(*,*) 'a'
00386         ELSE
00387 !           write(*,*) 'b', 1.2/(1.2**i)
00388            logtemp=logtemp-.9d0**i
00389         END IF
00390 !        IF (MPI_ID == 0) write(*,*) 'adjusted temp to', 10**logTemp, rho, IICoolingRate(10d0**logtemp), logtemp
00391      END DO
00392      GetIICoolEqTemp=10d0**logTemp
00393   END FUNCTION GetIICoolEqTemp
00394 
00395 
00397   FUNCTION GetIICoolEqDensity(p)
00398      REAL(KIND=qPrec) :: logrho,p, GetIICoolEqDensity
00399      INTEGER :: i
00400      logrho=7d0
00401      !using .1 allows for a range of +- 10 in log space
00402      !If we want 16 digits of accuracy we need -16/log10(.9)=360
00403       DO i=1, 350 !raise lower bounds to avoid missing local minima
00404         IF (10d0**logrho*IIHeatingRate() - (10d0**logrho)**2*IICoolingRate(p/10d0**logrho)  > 0d0) THEN
00405            logrho=logrho-.9d0**i
00406 !           write(*,*) 'a', 10d0**logrho,10d0**logrho*Y * ( 1d0 - 10d0**logrho*IICoolingRate(p/10d0**logrho))
00407         ELSE
00408            logrho=logrho+.9d0**i
00409 !           write(*,*) 'b', 10d0**logrho,10d0**logrho*Y * ( 1d0 - 10d0**logrho*IICoolingRate(p/10d0**logrho))
00410         END IF
00411 !        IF (MPI_ID == 0) write(*,*) 'adjusted temp to', 10**logTemp, rho, IICoolingRate(10d0**logtemp), logtemp
00412      END DO
00413      GetIICoolEqDensity=10d0**logrho
00414    END FUNCTION GetIICoolEqDensity
00415 
00416   FUNCTION IICoolingRate(T)
00417      REAL(KIND=qPrec) :: T,IICoolingRate
00418      IICoolingRate=IICoolPar(1)*(IICoolPar(2)*exp(-IICoolPar(3)/(T+IICoolPar(4)))+IICoolPar(5)*sqrt(T)*exp(-IICoolPar(6)/max(1d0,T+IICoolPar(7))))
00419   END FUNCTION IICoolingRate
00420 
00421 
00422   FUNCTION IIHeatingRate()
00423      REAL(KIND=qPrec) :: IIHeatingRate
00424      IIHeatingRate=IICoolPar(1)
00425    END FUNCTION IIHeatingRate
00426 
00433   SUBROUTINE NEQCooling(q,dqdt,coolingobj,lform)
00434     REAL(KIND=qPREC) :: q(:), dqdt(:)
00435     TYPE(CoolingDef), POINTER :: coolingobj
00436     LOGICAL :: lform
00437 
00438     IF (lform .eqv. PRIMITIVE) THEN
00439        dqdt(iE) = dqdt(iE) - NEQCoolingStrength(q,coolingobj,lform) * gamma1
00440     ELSE
00441        dqdt(iE) = dqdt(iE) - NEQCoolingStrength(q,coolingobj,lform)
00442     END IF
00443 
00444     dqdt(nSpeciesLo:nSpeciesHi) = dqdt(nSpeciesLo:nSpeciesHi) + NEQSpeciesChange(q,lform) 
00445 
00446   END SUBROUTINE NEQCooling
00447 
00448   FUNCTION NEQCoolingStrength(q,coolingobj,lform)
00449      REAL(KIND=qPrec) :: q(:), NEQCoolingStrength
00450      TYPE(CoolingDef), POINTER :: coolingobj
00451      LOGICAL :: lform
00452      REAL(KIND=qPrec) :: P, T0, mu
00453      REAL(KIND=qPrec), DIMENSION(0:nSpeciesHi) :: nvec
00454      REAL(KIND=qPrec), DIMENSION(NrHydroVars) :: f
00455      f = 0d0
00456 
00457      P = MERGE(PrimPress(q),Press(q),lform .eqv. PRIMITIVE)
00458      CALL GetNEQvars(q, mu, nvec)               ! defined in neqcooling.f90
00459      T0 = P/q(1) * TempScale * mu
00460   
00461      CALL Cool_Derivatives(q,f,T0,nvec)         ! defined in i_evolve.f90
00462 
00463      IF (T0 < coolingobj%floortemp) THEN
00464         NEQCoolingStrength = 0d0
00465      ELSE
00466         NEQCoolingStrength = f(iE)              ! scaling is done in i_evolve.f90
00467      END IF
00468 
00469   END FUNCTION NEQCoolingStrength
00470 
00471 
00472 
00480   SUBROUTINE ZCooling(q,dqdt,coolingobj,lform)
00481     REAL(KIND=qPREC) :: q(:), dqdt(:)
00482     TYPE(CoolingDef), POINTER :: coolingobj
00483     LOGICAL :: lform
00484 
00485     IF (lform .eqv. PRIMITIVE) THEN
00486        dqdt(iE) = dqdt(iE) - ZCoolingStrength(q,coolingobj,lform) * gamma1
00487     ELSE
00488        dqdt(iE) = dqdt(iE) - ZCoolingStrength(q,coolingobj,lform)
00489     END IF
00490 
00491     dqdt(nSpeciesLo:nSpeciesHi) = dqdt(nSpeciesLo:nSpeciesHi) + NEQSpeciesChange(q,lform) 
00492 
00493   END SUBROUTINE ZCooling
00494 
00495   FUNCTION ZCoolingStrength(q,coolingobj,lform)
00496      REAL(KIND=qPrec) :: q(:), ZCoolingStrength
00497      TYPE(CoolingDef),POINTER :: coolingobj
00498      REAL(KIND=qPrec) :: P,T0,nH,ne,x,mu
00499      LOGICAL :: lform
00500      REAL(KIND=qPrec), DIMENSION(0:nSpeciesHi) :: nvec
00501      REAL(KIND=qPrec), DIMENSION(NrHydroVars) :: f
00502      f = 0d0
00503 
00504      P = MERGE(PrimPress(q),Press(q),lform .eqv. PRIMITIVE)
00505      CALL GetNEQvars(q, mu, nvec)               ! defined in neqcooling.f90
00506      T0 = P/q(1) * TempScale * mu
00507 
00508      CALL Cool_Derivatives(q,f,T0,nvec,ZCool)         ! defined in i_evolve.f90
00509  
00510      CALL GetZvars(nvec, ne, x, nH)
00511 
00512      IF (T0 < coolingobj%floortemp) THEN
00513         ZCoolingStrength = 0d0
00514      ELSE
00515         IF(ZCoolingRate(ne,T0,x) == 0d0) THEN
00516            ZCoolingStrength = f(iE)
00517         ELSE
00518            ZCoolingStrength = f(iE) + (1d0-Zweight(T0))*metal_loss + Zweight(T0) * nH**2d0 * ZCoolingRate(ne,T0,x) * coolingobj%ScaleCool
00519         END IF
00520      END IF
00521 
00522   END FUNCTION ZCoolingStrength
00523 
00524   ! weight function used to average metal_loss and ZCoolingRate so that ZCoolingRate is heavily weighted
00525   ! in the middle of the Z cooling table, and goes to zero at the table limits. metal_loss is from NEQ
00526   ! cooling and is weighted in the opposited direction of ZCoolingRate.
00527   FUNCTION Zweight(T)
00528      REAL(KIND=qPrec) :: Zweight, T, midT
00529            
00530      midT = tempmax - (tempmax+tempmin)/2d0
00531      Zweight = (1d0/midT)**2d0 * (tempmax-T)*(T-tempmin) 
00532 !     Zweight = 1d0   
00533  
00534   END FUNCTION Zweight
00535 
00536   FUNCTION NEQSpeciesChange(q,lform)
00537      REAL(KIND=qPrec) :: q(:), NEQSpeciesChange(nSpeciesLo:nSpeciesHi)
00538      LOGICAL :: lform
00539      REAL(KIND=qPrec) :: P,T0,mu
00540      REAL(KIND=qPrec), DIMENSION(0:nSpeciesHi) :: nvec
00541      REAL(KIND=qPrec), DIMENSION(NrHydroVars) :: f
00542      f = 0d0
00543 
00544      P = MERGE(PrimPress(q),Press(q),lform .eqv. PRIMITIVE)
00545      CALL GetNEQvars(q, mu, nvec)               ! defined in neqcooling.f90
00546      T0 = P/q(1) * TempScale * mu
00547 
00548      CALL Cool_Derivatives(q,f,T0,nvec)         ! defined in i_evolve.f90
00549 
00550      NEQSpeciesChange = f(nSpeciesLo:nSpeciesHi)
00551 
00552   END FUNCTION NEQSpeciesChange
00553 
00554   ! Gets electron density, ionization fraction, and hydrogen number density
00555   SUBROUTINE GetZvars(nvecin,neout,xout,nHout)
00556      REAL(KIND=qPrec), DIMENSION(0:nSpeciesHi), INTENT(IN) :: nvecin
00557      REAL(KIND=qPrec) :: nneuc, npart
00558      REAL(KIND=qPrec), INTENT(OUT) :: neout, xout, nHout
00559 
00560      nneuc = nvecin(iH2)+nvecin(iH)+nvecin(iHII)+nvecin(iHe)+nvecin(iHeII)+nvecin(iHeIII)
00561      neout = nvecin(iHII) + nvecin(iHeII) + 2d0*nvecin(iHeIII) + nmine*nneuc
00562      npart = nneuc + neout
00563      nHout = nvecin(iH) + nvecin(iHII)
00564      xout = nvecin(iHII)/nHout
00565 
00566   END SUBROUTINE GetZvars
00567 
00568   ! trilinear interpolation algorithm based on formula on Paul Bourke's website
00569   FUNCTION ZCoolingRate(nein,tempin,xin)
00570      REAL(KIND=qPrec) :: lognein,logxin,ZCoolingRate,rate1,rate2,rate3,rate4,rate5,rate6,rate7,rate8,lintempin
00571      REAL(KIND=qPrec), INTENT(IN) :: nein, tempin, xin
00572      REAL(KIND=qPrec) :: lognestep,tempstep,logxstep
00573      INTEGER :: lognelow, lognehigh, templow, temphigh, logxlow, logxhigh 
00574 
00575      lognein = LOG10(nein)
00576      logxin = LOG10(xin)
00577      lintempin = tempin
00578 
00579      ! Compute step sizes of table Zcooling.tab
00580      lognestep = (lognemax-lognemin)/REAL(nDensities-1,qPrec)
00581      tempstep = (tempmax-tempmin)/REAL(nTemps-1,qPrec)
00582      logxstep = (logxmax-logxmin)/REAL(nXs-1,qPrec)
00583 
00584      ! Scale "in" values to fit within integer ranges of ZCoolingTab
00585      lognein = (lognein-lognemin)/lognestep + 1d0
00586      lintempin = (lintempin-tempmin)/tempstep + 1d0
00587      logxin = (logxin-logxmin)/logxstep + 1d0
00588 
00589      ! If the point is outside table range, return 0 cooling (essentially reverts to only being NEQ cooling)
00590      IF(lognein < 1 .OR. lognein > nDensities .OR. lintempin < 1 .OR. lintempin > nTemps .OR. logxin < 1 .OR. logxin > nXs) THEN
00591         ZCoolingRate = 0d0
00592      ELSE
00593 
00594         ! Find integers nearest to (lognein,tempin,logxin)
00595         lognelow = FLOOR(lognein); lognehigh = CEILING(lognein);
00596         templow = FLOOR(lintempin); temphigh = CEILING(lintempin);
00597         logxlow = FLOOR(logxin); logxhigh = CEILING(logxin);
00598 
00599         ! Match 8 closest points with their corresponding rates
00600         rate1 = ZCoolingTab(lognelow,templow,logxlow)
00601         rate2 = ZCoolingTab(lognehigh,templow,logxlow)
00602         rate3 = ZCoolingTab(lognelow,temphigh,logxlow)
00603         rate4 = ZCoolingTab(lognelow,templow,logxhigh)
00604         rate5 = ZCoolingTab(lognehigh,templow,logxhigh)
00605         rate6 = ZCoolingTab(lognelow,temphigh,logxhigh)
00606         rate7 = ZCoolingTab(lognehigh,temphigh,logxlow)
00607         rate8 = ZCoolingTab(lognehigh,temphigh,logxhigh)
00608 
00609         ! With low values --> 0 and high values --> 1, scale "in" values accordingly
00610         lognein = lognein-lognelow; lintempin = lintempin-templow; logxin = logxin-logxlow;
00611 
00612         ! With scaled values, the 8 closest points become the vertices of a unit cube
00613         ! aligned with the origin, and the interpolation is much simpler
00614         ZCoolingRate = rate1*(1d0-lognein)*(1d0-lintempin)*(1d0-logxin)
00615         ZCoolingRate = ZCoolingRate + rate2*lognein*(1d0-lintempin)*(1d0-logxin)
00616         ZCoolingRate = ZCoolingRate + rate3*(1d0-lognein)*lintempin*(1d0-logxin)
00617         ZCoolingRate = ZCoolingRate + rate4*(1d0-lognein)*(1d0-lintempin)*logxin
00618         ZCoolingRate = ZCoolingRate + rate5*lognein*(1d0-lintempin)*logxin
00619         ZCoolingRate = ZCoolingRate + rate6*(1d0-lognein)*lintempin*logxin
00620         ZCoolingRate = ZCoolingRate + rate7*lognein*lintempin*(1d0-logxin)
00621         ZCoolingRate = ZCoolingRate + rate8*lognein*lintempin*logxin
00622 
00623      END IF
00624   END FUNCTION ZCoolingRate
00625 
00626   ! ==========================================
00627   ! =      Cooling creation/destruction      =
00628   ! =      and list manipulation functions   =
00629   ! ==========================================
00630 
00631   SUBROUTINE CreateCoolingObject(coolingobj,userid)
00632      ! Interface declarations
00633      !INTEGER :: dummy
00634      TYPE(CoolingDef),POINTER :: coolingobj
00635      INTEGER,OPTIONAL :: userid
00636 
00637 !     IF(ASSOCIATED(coolingobj)) THEN
00638 !        PRINT*,'cooling_source.f90::CreateCoolingObject error -- object already associated. Halting.'
00639 !        STOP
00640 !     END IF
00641 
00642      ALLOCATE(coolingobj)
00643      NULLIFY(coolingobj%previous)
00644      NULLIFY(coolingobj%next)
00645      IF(PRESENT(userid)) THEN
00646         coolingobj%id=userid
00647      ELSE
00648         iCoolingObjID=iCoolingObjID+1
00649         coolingobj%id=iCoolingObjID
00650      END IF
00651      coolingobj%lhaspriority=.true.
00652 
00653      ! setup defaults: analytic cooling T^0.5, conservative w/o diagnostic field
00654      coolingobj%lprimitive=.false.
00655      coolingobj%icooling=AnalyticCool
00656      coolingobj%alpha=1d0
00657      coolingobj%beta=5d-1
00658      coolingobj%floortemp=1d-1
00659      coolingobj%mintemp=1d-2
00660 
00661      CALL AddCoolingObjToList(coolingobj)
00662   END SUBROUTINE CreateCoolingObject
00663 
00664   SUBROUTINE DestroyCoolingObject(coolingobj,id)
00665      TYPE(CoolingDef),POINTER :: coolingobj
00666      INTEGER,OPTIONAL :: id
00667 
00668      IF(PRESENT(id)) THEN
00669         coolingobj=>firstcoolingobj
00670         DO WHILE(ASSOCIATED(coolingobj))
00671            IF(coolingobj%id==id) THEN
00672               EXIT
00673            ELSE
00674               coolingobj=>coolingobj%next
00675            END IF
00676         END DO
00677      END IF
00678      CALL RemoveCoolingObjFromList(coolingobj)
00679      DEALLOCATE(coolingobj)
00680      NULLIFY(coolingobj)
00681   END SUBROUTINE DestroyCoolingObject
00682 
00683   SUBROUTINE AddCoolingObjToList(coolingobj)
00684      TYPE(CoolingDef),POINTER :: coolingobj
00685 
00686      IF(.NOT. ASSOCIATED(firstcoolingobj)) THEN ! first cooling object only
00687         firstcoolingobj=>coolingobj
00688         lastcoolingobj=>coolingobj
00689      ELSE
00690         coolingobj%previous=>lastcoolingobj
00691         lastcoolingobj%next=>coolingobj
00692         lastcoolingobj=>coolingobj
00693      END IF
00694   END SUBROUTINE AddCoolingObjToList
00695 
00696   SUBROUTINE RemoveCoolingObjFromList(coolingobj)
00697      TYPE(CoolingDef),POINTER :: coolingobj
00698 
00699      IF(ASSOCIATED(coolingobj%previous)) THEN
00700         coolingobj%previous%next=>coolingobj%next
00701      ELSE
00702         firstcoolingobj=>coolingobj%next
00703      END IF
00704 
00705      IF(ASSOCIATED(coolingobj%next)) THEN
00706         coolingobj%next%previous=>coolingobj%previous
00707      ELSE
00708         lastcoolingobj=>coolingobj%previous
00709         NULLIFY(lastcoolingobj%next)
00710      END IF
00711   END SUBROUTINE RemoveCoolingObjFromList
00712 
00713   LOGICAL FUNCTION CoolingCheck()
00714      TYPE(CoolingDef),POINTER :: coolingobj
00715      coolingcheck=.FALSE.
00716      coolingobj=>firstcoolingobj
00717      DO WHILE(ASSOCIATED(coolingobj))
00718         ! We need just one of the objects to exist 
00719         ! to keep lsrc true
00720         IF(coolingobj%iCooling/=0) THEN
00721            !          IF(testregion(coolingobj,x,dx) .AND. testregime(coolingobj)) THEN
00722            CoolingCheck=.TRUE.
00723            !          END IF
00724         END IF
00725         coolingobj=>coolingobj%next
00726      END DO
00727   END FUNCTION CoolingCheck
00728 
00729 
00730   ! ==========================================
00731   ! =       Initialization functions         =
00732   ! ==========================================
00733 
00734   SUBROUTINE InitCoolingTracers(sources)
00735      ! Interface declarations
00736       TYPE(SourcesDef),POINTER :: sources
00737      ! Internal declarations
00738      TYPE(CoolingDef),POINTER :: coolingObj
00739      INTEGER :: i
00740      coolingobj=>firstcoolingobj
00741      DO WHILE(ASSOCIATED(coolingObj))
00742         SELECT CASE(coolingObj%iCooling)
00743         CASE(NoCool)
00744         CASE(AnalyticCool)
00745         CASE(DMCool)
00746         CASE(NEQCool)
00747            CALL Addneqtracers
00748         CASE(ZCool)
00749            CALL Addneqtracers
00750         END SELECT
00751         coolingobj=>coolingobj%next
00752      END DO
00753   END SUBROUTINE InitCoolingTracers
00754 
00755   SUBROUTINE InitCoolingElliptics(sources)
00756      ! Interface declarations
00757      TYPE(SourcesDef),POINTER :: sources
00758      ! Internal declarations
00759      TYPE(CoolingDef),POINTER :: coolingObj
00760      INTEGER :: i
00761      coolingobj=>firstcoolingobj
00762      DO WHILE(ASSOCIATED(coolingObj))
00763         SELECT CASE(coolingObj%iCooling)
00764         CASE(NoCool)
00765         CASE(AnalyticCool)
00766         CASE(DMCool)
00767         CASE(NEQCool)
00768         CASE(ZCool)
00769         END SELECT
00770         coolingobj=>coolingobj%next
00771      END DO
00772   END SUBROUTINE InitCoolingElliptics
00773 
00776   SUBROUTINE CoolingFinalizeInit
00777      TYPE(CoolingDef),POINTER :: coolingobj
00778 
00779      coolingobj=>firstcoolingobj
00780      DO WHILE(ASSOCIATED(coolingobj))
00781         SELECT CASE(coolingobj%iCooling)
00782         CASE(NoCool)
00783         CASE(AnalyticCool)
00784            CALL InitAnalyticCool(coolingobj)
00785         CASE(DMCool)
00786            CALL InitDMCool
00787            coolingobj%scalecool=nScale**2*TimeScale/pScale
00788         CASE(IICool)
00789            CALL InitIICool(coolingobj)
00790         CASE(NEQCool)
00791            CALL InitNeqCool
00792         CASE(ZCool)
00793            coolingobj%scalecool=TimeScale/pScale
00794            CALL InitNeqCool
00795            CALL InitZCool
00796         END SELECT
00797         coolingobj=>coolingobj%next
00798      END DO
00799   END SUBROUTINE CoolingFinalizeInit
00800 
00802   SUBROUTINE InitAnalyticCool(coolingobj)
00803      TYPE(CoolingDef),POINTER :: coolingobj
00804      !For power laws that have the same strength as brehm strahlung at a given shock speed (cs*Mach) 
00805      !use alpha=4.76e-20 (ergs*cm^3/s/K^.5) * (3d0/16d0*TempScale*(cs*Mach)**2)**(.5d0)  * (3d0/16d0 * (cs*Mach)**2) **(-coolingobj%beta)
00806 
00807      !Note if beta=.5 this is just pure Brehmstrahlung and
00808      ! alpha = 4.76e-20 (ergs*cm^3/s/K^.5) * (TempScale)**(.5d0)
00809      !See Imamura, Wolff, and Durisen (1983) although their alpha is our beta and their lambda is our alpha
00810 
00811      !alpha should have units of ergs*cm^3/s/K^beta
00812      coolingobj%scalecool=nScale**2/(pScale/Timescale)
00813      write(*,*) 'scalecool=', coolingobj%scalecool
00814   END SUBROUTINE InitAnalyticCool
00815 
00818   SUBROUTINE InitIICool(coolingobj)
00819      TYPE(CoolingDef),POINTER :: coolingobj
00820      !
00821      INTEGER :: i,nInput
00822      REAL :: dummy
00823 
00824      coolingobj%ScaleCool = 1d0/(pScale/TimeScale)
00825 
00826 !     Switched to analytic expression instead of using a table
00827 !     IF(ALLOCATED(IICoolingTab)) RETURN
00828 !     OPEN(UNIT=85,FILE='TABLES/IIcooling.tab')
00829 !     READ(85,*) nInput
00830 !     ALLOCATE(IICoolingTab(nInput))
00831 !     DO i=1,nInput
00832 !        READ(85,*) dummy, IICoolingTab(i)
00833 !     END DO
00834 !     CLOSE(85)
00835   END SUBROUTINE InitIICool
00836 
00843   SUBROUTINE InitZCool
00844 
00845      IF(ALLOCATED(ZCoolingTab)) RETURN
00846      OPEN(UNIT=85,FILE='TABLES/Zcooling.tab')
00847      READ(85,*) nDensities,nTemps,nXs,lognemin,lognemax,tempmin,tempmax,logxmin,logxmax
00848      ALLOCATE(ZCoolingTab(nDensities,nTemps,nXs))
00849      READ(85,*) ZCoolingTab
00850      CLOSE(85)
00851   END SUBROUTINE InitZCool
00852 
00856   SUBROUTINE Cooling_ReadObjectFromChombo(chandle, cooling_object)
00857 
00858      USE ChomboDeclarations, ONLY: ChomboHandle
00859      USE HDF5Declarations, ONLY: Read_Slab_From_Dataset_Int, Read_Slab_From_Dataset_Double
00860 
00861      TYPE(ChomboHandle), POINTER :: chandle
00862      TYPE(CoolingDef), POINTER :: cooling_object
00863 
00864      INTEGER, DIMENSION(1), TARGET :: int_buffer_array
00865      REAL(KIND=qPrec), DIMENSION(1), TARGET :: dbl_buffer_array
00866      INTEGER, DIMENSION(:), POINTER :: int_buffer
00867      REAL(KIND=qPrec), DIMENSION(:), POINTER :: dbl_buffer
00868 
00869 
00870      int_buffer => int_buffer_array
00871      dbl_buffer => dbl_buffer_array
00872 
00873      int_buffer = 0
00874      dbl_buffer = 0.d0
00875 
00876      CALL Read_Slab_From_Dataset_Int("iCooling", &
00877           chandle%source_group_id, &
00878           int_buffer, &
00879           chandle%source_offset)
00880 
00881      cooling_object%iCooling = int_buffer(1)
00882 
00883      CALL Read_Slab_From_Dataset_Double("alpha", &
00884           chandle%source_group_id, &
00885           dbl_buffer, &
00886           chandle%source_offset)
00887 
00888      cooling_object%alpha = dbl_buffer(1)
00889 
00890      CALL Read_Slab_From_Dataset_Double("beta", &
00891           chandle%source_group_id, &
00892           dbl_buffer, &
00893           chandle%source_offset)
00894 
00895      cooling_object%beta = dbl_buffer(1)
00896 
00897      CALL Read_Slab_From_Dataset_Int("lPrimitive", &
00898           chandle%source_group_id, &
00899           int_buffer, &
00900           chandle%source_offset)
00901 
00902      cooling_object%lPrimitive = (int_buffer(1) /= 0)
00903 
00904      ! Write position components to datasets.  Be sure to index the arrays using the i:i 
00905      ! notation, as this passes a 1-element array to the subroutine instead of a scalar.
00906      CALL Read_Slab_From_Dataset_Double("cartesian_xlower", &
00907           chandle%source_group_id, &
00908           dbl_buffer, &
00909           chandle%source_offset)
00910 
00911      cooling_object%xlower(1) = dbl_buffer(1)
00912 
00913      CALL Read_Slab_From_Dataset_Double("cartesian_ylower", &
00914           chandle%source_group_id, &
00915           dbl_buffer, &
00916           chandle%source_offset)
00917 
00918      cooling_object%xlower(2) = dbl_buffer(1)
00919 
00920      CALL Read_Slab_From_Dataset_Double("cartesian_zlower", &
00921           chandle%source_group_id, &
00922           dbl_buffer, &
00923           chandle%source_offset)
00924 
00925      cooling_object%xlower(3) = dbl_buffer(1)
00926 
00927      CALL Read_Slab_From_Dataset_Double("cartesian_xupper", &
00928           chandle%source_group_id, &
00929           dbl_buffer, &
00930           chandle%source_offset)
00931 
00932      cooling_object%xupper(1) = dbl_buffer(1)
00933 
00934      CALL Read_Slab_From_Dataset_Double("cartesian_yupper", &
00935           chandle%source_group_id, &
00936           dbl_buffer, &
00937           chandle%source_offset)
00938 
00939      cooling_object%xupper(2) = dbl_buffer(1)
00940 
00941      CALL Read_Slab_From_Dataset_Double("cartesian_zupper", &
00942           chandle%source_group_id, &
00943           dbl_buffer, &
00944           chandle%source_offset)
00945 
00946      cooling_object%xupper(3) = dbl_buffer(1)
00947 
00948      CALL Read_Slab_From_Dataset_Double("spherical_pos_x", &
00949           chandle%source_group_id, &
00950           dbl_buffer, &
00951           chandle%source_offset)
00952 
00953      cooling_object%pos(1) = dbl_buffer(1)
00954 
00955      CALL Read_Slab_From_Dataset_Double("spherical_pos_y", &
00956           chandle%source_group_id, &
00957           dbl_buffer, &
00958           chandle%source_offset)
00959 
00960      cooling_object%pos(3) = dbl_buffer(1)
00961 
00962      CALL Read_Slab_From_Dataset_Double("spherical_pos_z", &
00963           chandle%source_group_id, &
00964           dbl_buffer, &
00965           chandle%source_offset)
00966 
00967      cooling_object%pos(3) = dbl_buffer(1)
00968 
00969      CALL Read_Slab_From_Dataset_Double("spherical_radius", &
00970           chandle%source_group_id, &
00971           dbl_buffer, &
00972           chandle%source_offset)
00973 
00974      cooling_object%radius = dbl_buffer(1)
00975 
00976      CALL Read_Slab_From_Dataset_Double("ScaleCool", &
00977           chandle%source_group_id, &
00978           dbl_buffer, &
00979           chandle%source_offset)
00980 
00981      cooling_object%ScaleCool = dbl_buffer(1)
00982 
00983      CALL Read_Slab_From_Dataset_Int("lHasPriority", &
00984           chandle%source_group_id, &
00985           int_buffer, &
00986           chandle%source_offset)
00987 
00988      cooling_object%lHasPriority = (int_buffer(1) /= 0)
00989 
00990      chandle%source_offset = chandle%source_offset + 1
00991 
00992   END SUBROUTINE Cooling_ReadObjectFromChombo
00993 
00994 
00998   SUBROUTINE Cooling_InitChomboDatasets(chandle, obj_count)
00999 
01000      USE ChomboDeclarations, ONLY: ChomboHandle
01001      USE HDF5Declarations, ONLY: Initialize_HDF5_Dataset_Int, Initialize_HDF5_Dataset_Double
01002 
01003      TYPE(ChomboHandle), POINTER :: chandle
01004      INTEGER :: obj_count
01005 
01006      TYPE(CoolingDef), POINTER :: cooling_object
01007      INTEGER :: i_err
01008      INTEGER :: iFixed
01009 
01010 
01011      IF (.NOT. ASSOCIATED(chandle)) THEN
01012         PRINT *, "Cooling_InitChomboDatasets error::invalid Chombo handle."
01013         STOP
01014      END IF
01015 
01016      IF (obj_count < 0) THEN
01017         PRINT *, "Cooling_InitChomboDatasets error::invalid object count ", obj_count, "."
01018         STOP
01019      END IF
01020 
01021      ! Position coordinates
01022      CALL Initialize_HDF5_Dataset_Int("iCooling", chandle%source_group_id, obj_count)
01023 
01024      CALL Initialize_HDF5_Dataset_Double("alpha", chandle%source_group_id, obj_count)
01025      CALL Initialize_HDF5_Dataset_Double("beta", chandle%source_group_id, obj_count)
01026 
01027      ! Logical variable.
01028      CALL Initialize_HDF5_Dataset_Int("lPrimitive", chandle%source_group_id, obj_count)
01029 
01030      CALL Initialize_HDF5_Dataset_Double("FloorTemp", chandle%source_group_id, obj_count)
01031      CALL Initialize_HDF5_Dataset_Double("MinTemp", chandle%source_group_id, obj_count)
01032 
01033      ! Cartesian region coordinates, stored in the same fashion as the sink particle coordinates.
01034      ! I don't know if this is helpful or not.
01035      CALL Initialize_HDF5_Dataset_Double("cartesian_xlower", chandle%source_group_id, obj_count)
01036      CALL Initialize_HDF5_Dataset_Double("cartesian_xupper", chandle%source_group_id, obj_count)
01037      CALL Initialize_HDF5_Dataset_Double("cartesian_ylower", chandle%source_group_id, obj_count)
01038      CALL Initialize_HDF5_Dataset_Double("cartesian_yupper", chandle%source_group_id, obj_count)
01039      CALL Initialize_HDF5_Dataset_Double("cartesian_zlower", chandle%source_group_id, obj_count)
01040      CALL Initialize_HDF5_Dataset_Double("cartesian_zupper", chandle%source_group_id, obj_count)
01041 
01042      ! Spherical region coordinates, stored in the same fashion as the sink particle coordinates.
01043      ! I don't know if this is helpful or not.
01044      CALL Initialize_HDF5_Dataset_Double("spherical_pos_x", chandle%source_group_id, obj_count)
01045      CALL Initialize_HDF5_Dataset_Double("spherical_pos_y", chandle%source_group_id, obj_count)
01046      CALL Initialize_HDF5_Dataset_Double("spherical_pos_z", chandle%source_group_id, obj_count)
01047 
01048      ! Radius of spherical region (scalar).
01049      CALL Initialize_HDF5_Dataset_Double("spherical_radius", chandle%source_group_id, obj_count)
01050 
01051      CALL Initialize_HDF5_Dataset_Double("ScaleCool", chandle%source_group_id, obj_count)
01052 
01053      ! Logical variables.
01054      CALL Initialize_HDF5_Dataset_Int("lHasPriority", chandle%source_group_id, obj_count)
01055 
01056   END SUBROUTINE Cooling_InitChomboDatasets
01057 
01060   SUBROUTINE Cooling_WriteObjectToChombo(chandle, cooling_object)
01061 
01062      USE ChomboDeclarations, ONLY: ChomboHandle
01063      USE HDF5Declarations, ONLY: Write_Slab_To_Dataset_Int, Write_Slab_To_Dataset_Double
01064 
01065      TYPE(ChomboHandle), POINTER :: chandle
01066      TYPE(CoolingDef), POINTER :: cooling_object
01067 
01068 
01069      CALL Write_Slab_To_Dataset_Int("iCooling", &
01070           chandle%source_group_id, &
01071           (/ cooling_object%iCooling /), &
01072           chandle%source_offset)
01073 
01074      CALL Write_Slab_To_Dataset_Double("alpha", &
01075           chandle%source_group_id, &
01076           (/ cooling_object%alpha /), &
01077           chandle%source_offset)
01078 
01079      CALL Write_Slab_To_Dataset_Double("beta", &
01080           chandle%source_group_id, &
01081           (/ cooling_object%beta /), &
01082           chandle%source_offset)
01083 
01084      CALL Write_Slab_To_Dataset_Int("lPrimitive", &
01085           chandle%source_group_id, &
01086           (/ BoolToInt(cooling_object%lPrimitive) /), &
01087           chandle%source_offset)
01088 
01089      ! Write position components to datasets.  Be sure to index the arrays using the i:i 
01090      ! notation, as this passes a 1-element array to the subroutine instead of a scalar.
01091      CALL Write_Slab_To_Dataset_Double("cartesian_xlower", &
01092           chandle%source_group_id, &
01093           cooling_object%xlower(1:1), &
01094           chandle%source_offset)
01095 
01096      CALL Write_Slab_To_Dataset_Double("cartesian_ylower", &
01097           chandle%source_group_id, &
01098           cooling_object%xlower(2:2), &
01099           chandle%source_offset)
01100 
01101      CALL Write_Slab_To_Dataset_Double("cartesian_zlower", &
01102           chandle%source_group_id, &
01103           cooling_object%xlower(3:3), &
01104           chandle%source_offset)
01105 
01106      CALL Write_Slab_To_Dataset_Double("cartesian_xupper", &
01107           chandle%source_group_id, &
01108           cooling_object%xupper(1:1), &
01109           chandle%source_offset)
01110 
01111      CALL Write_Slab_To_Dataset_Double("cartesian_yupper", &
01112           chandle%source_group_id, &
01113           cooling_object%xupper(2:2), &
01114           chandle%source_offset)
01115 
01116      CALL Write_Slab_To_Dataset_Double("cartesian_zupper", &
01117           chandle%source_group_id, &
01118           cooling_object%xupper(3:3), &
01119           chandle%source_offset)
01120 
01121      CALL Write_Slab_To_Dataset_Double("spherical_pos_x", &
01122           chandle%source_group_id, &
01123           cooling_object%pos(1:1), &
01124           chandle%source_offset)
01125 
01126      CALL Write_Slab_To_Dataset_Double("spherical_pos_y", &
01127           chandle%source_group_id, &
01128           cooling_object%pos(2:2), &
01129           chandle%source_offset)
01130 
01131      CALL Write_Slab_To_Dataset_Double("spherical_pos_z", &
01132           chandle%source_group_id, &
01133           cooling_object%pos(3:3), &
01134           chandle%source_offset)
01135 
01136      CALL Write_Slab_To_Dataset_Double("spherical_radius", &
01137           chandle%source_group_id, &
01138           (/ cooling_object%radius /), &
01139           chandle%source_offset)
01140 
01141      CALL Write_Slab_To_Dataset_Double("ScaleCool", &
01142           chandle%source_group_id, &
01143           (/ cooling_object%ScaleCool /), &
01144           chandle%source_offset)
01145 
01146      CALL Write_Slab_To_Dataset_Int("lHasPriority", &
01147           chandle%source_group_id, &
01148           (/ BoolToInt(cooling_object%lHasPriority) /), &
01149           chandle%source_offset)
01150 
01151      chandle%source_offset = chandle%source_offset + 1
01152 
01153   END SUBROUTINE Cooling_WriteObjectToChombo
01154 END MODULE CoolingSrc
 All Classes Files Functions Variables