Scrambler  1
BonnorEbertSphere/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 BonnorEbertSphere 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 !#########################################################################
00023 !============================================================================
00024 ! This problem places a bonner ebert sphere in the center of the grid. 
00025 !============================================================================
00026 
00027 
00028 MODULE PROBLEM
00029   USE DataDeclarations
00030   USE GlobalDeclarations
00031   USE PhysicsDeclarations
00032   USE BE_MODULE
00033   USE WINDS
00034   USE CLUMPS
00035   USE Ambients
00036   USE Totals
00037 
00038   IMPLICIT NONE    ! It's safer to require explicit declarations
00039   SAVE             ! Save module information
00040   PRIVATE
00041 
00042   PUBLIC ProblemModuleInit, ProblemGridInit, &
00043        ProblemBeforeStep, ProblemAfterStep, ProblemSetErrFlag, ProblemBeforeGlobalStep
00044 
00045   TYPE(WindDef), POINTER :: Wind
00046   TYPE(ClumpDef), POINTER :: MyClump
00047   TYPE(AmbientDef), POINTER :: Ambient
00048   TYPE(TotalDef), POINTER :: Total
00049   INTEGER :: nWinds
00050   REAL(KIND=qPrec), DIMENSION(3) :: xloc=(/0,0,0/) !Clump Location -- Erica 1/18/2012
00051   REAL(KIND=qPrec) :: central_rho, clump_rad, rho_out, m_star, at, rho_c, gravo_thermal, time_f, iso_t, rho_weight! ,p_weight, tff, omega
00052   REAL :: xi
00053   REAL(KIND=qPrec) :: rhoOut, pOut, vxOut, vyOut, vzOut, BxOut, ByOut, BzOut, crit_rad, p_crit
00054   REAL(KIND=qPrec) :: Omega, tff
00055   
00056   NAMELIST /ProblemData/ xloc, central_rho, xi, clump_rad, rho_weight, omega, CellsPerJeansLength !p_weight -- var not used anymore. Also, remember when define omega in prob.data to * with timescale
00057   NAMELIST /AmbientData/ rhoOut, pOut, vxOut, vyOut, vzOut, BxOut, ByOut, BzOut
00058 
00059 
00060 CONTAINS
00061 
00062   SUBROUTINE ProblemModuleInit 
00063 
00064   INTEGER :: i, edge
00065 
00066     ! Read in problem params and vars
00067 
00068     OPEN(UNIT = PROBLEM_DATA_HANDLE, FILE = 'problem.data')
00069     READ(PROBLEM_DATA_HANDLE, NML=ProblemData)
00070     READ(PROBLEM_DATA_HANDLE, NML=AmbientData)
00071     CLOSE(PROBLEM_DATA_HANDLE)
00072     central_rho = central_rho*hmass*Xmu ! now is g/cc
00073     CALL CalcAmbientParams(central_rho, xi, clump_rad, M_star, at, rho_out, iso_t)
00074     CALL CreateAmbient(Ambient)
00075     !gravo_thermal = Sqrt((4*pi*G*central_rho)/(at**2))    ! Gravo_thermal radius -- constant that scales r to xi
00076     pOut=(rho_out/rscale)*(iso_t/tempscale)!((10*(aT**8))/((G**3)*(M_star**2)))/pscale
00077     rhoOut=(rho_out/rscale)/rho_weight
00078     Ambient%density=rhoOut
00079     Ambient%pressure=pOut
00080     Ambient%B(:)=(/BxOut,ByOut,BzOut/)
00081     Ambient%velocity(:)=(/vxOut, vyOut, vzOut/)
00082     crit_rad = (0.41*G*M_star)/at**2
00083     p_crit = ((1.4*(aT**8))/((G**3)*(M_star**2)))/pscale
00084     CALL CreateClump(MyClump) 
00085     MyClump%position = xloc
00086     MyClump%density = central_rho/rscale
00087     MyClump%density_profile = BE_PROFILE
00088     MyClump%temperature = iso_t/tempScale
00089     MyClump%omega = omega*72046714255133.6 !* timeScale to put omega in cu
00090     MyClump%radius = clump_rad/lscale
00091     CALL UpdateClump(MyClump)
00092       
00093       CALL CreateTotal(Total)
00094       Total%Field%Component=GASCOMP
00095       Total%Field%id=ivx
00096       Total%Field%name='Gas Px'
00097       
00098       CALL CreateTotal(Total)
00099       Total%Field%Component=PARTICLECOMP
00100       Total%Field%id=ivx
00101       Total%Field%name='Particle Px'
00102       
00103       CALL CreateTotal(Total)
00104       Total%Field%Component=BOTHCOMP
00105       Total%Field%id=ivx
00106       Total%Field%name='Combined Px'
00107 
00108       CALL CreateTotal(Total)
00109       Total%Field%Component=GASCOMP
00110       Total%Field%id=1
00111       Total%Field%name='Gas rho'
00112       
00113       CALL CreateTotal(Total)
00114       Total%Field%Component=PARTICLECOMP
00115       Total%Field%id=1
00116       Total%Field%name='Particle mass'
00117       
00118       CALL CreateTotal(Total)
00119       Total%Field%Component=BOTHCOMP
00120       Total%Field%id=1
00121       Total%Field%name='Combined mass'   
00122 
00123 
00124     tff= (0.5427)/(Sqrt(G*central_rho))
00125 
00126     IF (MPI_ID == 0) THEN 
00127        write(*,*) ambient%density, central_rho*BE_RHO(xi), rho_out, xi
00128        WRITE (*,*) 'The mass (in solar masses) of your star is:' , M_star/msolar 
00129        WRITE (*,*) 'The isothermal sound speed (cm/s) of your star and surrounding medium is:' , aT
00130        WRITE (*,*) 'The isothermal temperature (K) of your star and surrounding medium is:', iso_T
00131 
00132        IF ((central_rho/rho_out) < 14.1) THEN  
00133           WRITE (*,*) 'Your BE sphere IS stable. Rho_c/Rho_o = ', central_rho/rho_out, "(less than 14.1)"
00134        ELSE 
00135           WRITE(*,*) '!!Your BE sphere is NOT stable!! Rho_c/Rho_o =', central_rho/rho_out, "(greater than 14.1)"
00136        END IF
00137        time_f = clump_rad/at !sound crossing time in cgs
00138        time_f = time_f/timeScale !sound crossing time in cu
00139        WRITE(*,*) 'The sound crossing time for your clump in computational units is:', time_f
00140        WRITE(*,*) 'The central density of your clump (in CGS) is:', central_rho
00141        WRITE(*,*) 'The outter density of your clump (in CGS) is:', rho_out
00142        WRITE(*,*) 'The free fall time in seconds is', tff
00143        WRITE(*,*) 'G', g, 'Scalegrav', scalegrav
00144        WRITE(*,*) 'Omega in cgs', (0.1)/tff
00145     END IF
00146 
00147       nWinds=0
00148       DO i=1,nDim
00149          DO edge=1,2
00150             IF (Gmthbc(i,edge) == 1) THEN 
00151                nWinds=nWinds+1
00152                CALL CreateWind(Wind)
00153                Wind%dir=i
00154                Wind%density=Ambient%density
00155                Wind%temperature=Ambient%pressure/Ambient%density
00156                Wind%edge=edge
00157             END IF
00158          END DO
00159       END DO
00160 
00161 
00162   END SUBROUTINE ProblemModuleInit
00163 
00164   SUBROUTINE ProblemGridInit(Info) !Specific for each info structure
00165     Type (InfoDef) :: Info
00166     INTEGER        :: i, j, k, mx, my, mz
00167     INTEGER        :: mbc, rmbc, zrmbc
00168     INTEGER        :: nx, ny, nz
00169     INTEGER        :: ndim
00170     REAL(KIND=xPrec):: x, y, z, xlower, ylower, zlower,dx, dy, dz, r
00171 
00172 
00173 !To initialize collapse, uncomment the below line to 
00174 !lower the internal energy everywhere in the grid.
00175 !Since objects are initialized before the grid is, placing 
00176 !this line here effectively moves the clump out of equilibrium
00177 
00178 
00179 !Info%q(:,:,:,1) = (1.1)*Info%q(:,:,:,1) !UNCOMMENT TO INDUCE COLLAPSE
00180 !Info%q(:,:,:,iE) = (1.1)*Info%q(:,:,:,iE) !UNCOMMENT TO INDUCE COLLAPSE
00181 
00182 
00183 
00184 
00185 RETURN  
00186   
00187     ! Assigning shorthand notation for variables already defined in Info  
00188     rmbc=levels(Info%level)%gmbc(levels(Info%level)%step)! Sets up ghost cells and boundaries
00189     mx=Info%mx(1);  xlower=Info%Xbounds(1,1) 
00190     my=Info%mx(2);  ylower=Info%Xbounds(2,1)
00191     mz=Info%mx(3);  zlower=Info%Xbounds(3,1)
00192 
00193     dx=levels(Info%level)%dx
00194 
00195     ! These seem to be "nested" variable declarations a=b=c...?
00196 
00197 
00198     SELECT CASE(nDim)
00199     CASE(2)
00200        zrmbc=0
00201        mz=1
00202        zlower=0
00203        dz=0  
00204     CASE(3)
00205        zrmbc=rmbc
00206     END SELECT
00207 
00208     ! Initializing the environment
00209     !Info%q(:,:,:,:)=0.d0                                  ! Sanity check -- setting grid to zero
00210     !gravo_thermal = Sqrt((4*pi*G*central_rho)/(at**2))    ! Gravo_thermal radius -- constant that scales r to xi
00211     ! Must I allocate an array anywhere for Astrobear 2.0?
00212     ! i,j,k are indices of cells, and x,y,z are spatial cartesian coords
00213     
00214     Do k=1-zrmbc, mz+zrmbc
00215        Do j=1-rmbc, my+rmbc
00216           Do i=1-rmbc, mx+rmbc
00217              x=(xlower + (REAL(i, xPrec)-half)*dx)
00218              y=(ylower + (REAL(j, xPrec)-half)*dx)
00219              z=(zlower + (REAL(k, xPrec)-half)*dx) 
00220              r = Sqrt(x**2+y**2+z**2) 
00221              r = (r*lscale)*gravo_thermal                 ! takes r to cgs to xi 
00222              IF (r <= xi) THEN                            ! if r<pc
00223                 Info%q(i,j,k,1) = (central_rho*BE_RHO(REAL(r)))/rscale 
00224              ELSE 
00225                 Info%q(i,j,k,1) = rho_out/rscale                                              
00226              END IF
00227              ! M_star = 8*PI*central_rho*clump_rad**3 *((1 + alpha)/(6 + xi**2) - (6*alpha)/(6 + xi**2)**2 - alpha/(3**(1/alpha)*12*e + xi**2) + D*(2**(-alpha)-1)/ (1+2**(-alpha) * D* xi**2)*(1+D*xi**2))
00228              Info%q(i,j,k,iE) = (((Info%q(i,j,k,1)*rscale*(aT**2))/gamma)*gamma7)/pScale
00229           END DO
00230        END DO
00231     END DO
00232 
00233 
00234   END SUBROUTINE ProblemGridInit
00235 
00236   ! Applies Boundary Conditions
00237   SUBROUTINE ProblemBeforeStep(INFO)
00238 
00239     TYPE(INFODEF) :: INFO
00240     INTEGER :: i
00241 
00242   END SUBROUTINE ProblemBeforeStep
00243 
00244   SUBROUTINE ProblemAfterStep(INFO)
00245     TYPE(INFODEF) :: INFO
00246     ! No special after step instructions needed; this is a stub.
00247 
00248   END SUBROUTINE ProblemAfterStep
00249 
00250   SUBROUTINE ProblemSetErrFlag(INFO)
00251     TYPE(INFODEF) :: INFO
00252     ! No special instructions needed; this is a stub.
00253 
00254   END SUBROUTINE ProblemSetErrFlag
00255 
00256   SUBROUTINE ProblemBeforeGlobalStep(n)
00257      INTEGER :: n
00258   END SUBROUTINE ProblemBeforeGlobalStep
00259 
00260 
00261 END MODULE PROBLEM
 All Classes Files Functions Variables