Scrambler
1
|
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 StreamDisk 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 USE DataDeclarations 00037 USE PhysicsDeclarations 00038 USE ProcessingDeclarations 00039 USE Ambients 00040 USE Disks 00041 USE PointGravitySrc 00042 USE Profiles 00043 USE CoolingSrc 00044 USE Refinements 00045 USE Projections 00046 IMPLICIT NONE 00047 SAVE 00048 00049 PUBLIC ProblemModuleInit, ProblemGridInit, ProblemBeforeStep, ProblemAfterStep, ProblemSetErrFlag, ProblemBeforeGlobalStep 00050 00051 PRIVATE 00052 REAL(KIND=qPREC) :: ambient_density, ambient_temperature, !ambient params 00053 disk_softlength, disk_density, disk_temperature, disk_radius, disk_smoothing, disk_height, !disk params 00054 stream_density, stream_temperature, stream_velocity, stream_radius, stream_halfwidth, stream_halfheight, stream_xpos, !stream params 00055 central_mass, !star params 00056 hotspot_width, hotspot_pos(3) !hot spot params 00057 INTEGER :: disk_maxlevel !maximum level to refine disk 00058 LOGICAL :: lCooling,lGlobalSim, lMatchPressure, lRefineHotSpot !Logical switches 00059 INTEGER :: stream_tracer=0 !module integer handle for stream tracer 00060 TYPE(AmbientDef), POINTER :: Ambient !module pointer for ambient 00061 00062 CONTAINS 00063 00065 SUBROUTINE ProblemModuleInit() 00066 ! Declare Objects initialized here 00067 TYPE(DiskDef), POINTER :: Disk 00068 TYPE(PointGravityDef), POINTER :: PointGravity 00069 TYPE(CoolingDef),POINTER :: coolingobj 00070 TYPE(RefinementDef), POINTER :: DiskRefinement, HotSpotRefinement, DiskAntiRefinement, StreamInjectionRefinement, StreamRefinement 00071 TYPE(ProjectionDef), POINTER :: Projection 00072 INTEGER :: i, j, k ! Miscellaneous counters 00073 INTEGER :: nPoints=400 !number of points in reconstructing hydrostatic equilibrium solution 00074 REAL(KIND=qPREC) :: r, r_max, p_inf, phi, rho 00075 00076 NAMELIST /ProblemData/ ambient_density, ambient_temperature, & 00077 disk_softlength, disk_density, disk_temperature, disk_radius, disk_height, disk_smoothing, disk_maxlevel, & 00078 stream_density, stream_temperature, stream_velocity, stream_halfwidth, stream_halfheight, stream_xpos, & 00079 central_mass, & 00080 lRefineHotSpot, hotspot_width, hotspot_pos, & 00081 lCooling, lGlobalSim, lMatchPressure 00082 00083 OPEN(UNIT=PROBLEM_DATA_HANDLE, FILE='problem.data', STATUS="OLD") 00084 READ(PROBLEM_DATA_HANDLE,NML=ProblemData) 00085 CLOSE(PROBLEM_DATA_HANDLE) 00086 00087 00088 !Point Gravity Stuff 00089 IF (.NOT. lRestart) THEN 00090 CALL CreatePointGravityObject(PointGravity) 00091 ELSE 00092 PointGravity=>FirstPointGravityObj 00093 END IF 00094 PointGravity%soft_length=disk_softlength*disk_radius/lScale 00095 PointGravity%soft_function=SPLINESOFT 00096 PointGravity%mass=central_mass*mSolar/mScale 00097 00098 GravityDim = 3 !We want point gravity objects to have 1/r^2 force - not 1/r 00099 00100 ! Check that softening length is resolved to avoid lots of numerical dissipation and heating at disk center 00101 IF (PointGravity%soft_length/levels(0)%dx < 4d0) THEN 00102 IF (MPI_ID == 0) WRITE(*,*) 'PointGravity Softening Radius is only ', PointGravity%soft_length/levels(0)%dx, ' root cells' 00103 STOP 00104 END IF 00105 00106 00107 !Ambient Stuff 00108 CALL CreateAmbient(Ambient) 00109 Ambient%density=ambient_density/rScale 00110 Ambient%pressure=Ambient%density*ambient_temperature/TempScale 00111 Ambient%PersistInBoundaries=(gmthbc == 1) !Any open boundary should maintain ambient values 00112 00113 !Disk Stuff 00114 CALL CreateDisk(Disk) 00115 Disk%density=disk_density/rScale 00116 Disk%pressure=Disk%density*disk_temperature/TempScale 00117 Disk%radius=disk_radius/lScale 00118 Disk%soft_length=PointGravity%soft_length 00119 Disk%soft_function=PointGravity%soft_function 00120 Disk%central_mass=PointGravity%mass 00121 Disk%thickness=disk_smoothing 00122 Disk%height=disk_height 00123 Disk%PersistInBoundaries(1,1)=.NOT. lGlobalSim 00124 CALL AddTracer(Disk%iTracer, 'Disk Tracer') 00125 CALL Updatedisk(Disk) 00126 00127 00128 !Now we need to calculate maximum distance from disk to any given cell 00129 !First get maximum distance to each corner from (0d0,0d0,0d0) 00130 r_max=0d0 00131 DO i=1,2 00132 DO j=1,merge(2,1,nDim >= 2) 00133 DO k=1,merge(2,1,nDim >= 3) 00134 r_max=max(r_max, sqrt(GxBounds(1,i)**2+GxBounds(2,j)**2+GxBounds(3,k)**2)) 00135 END DO 00136 END DO 00137 END DO 00138 !Add in extra ghost zones 00139 r_max=r_max+sqrt(REAL(nDim))*levels(0)%dx*6 !make sure we have enough distance to include ghost zones... 00140 IF (lMatchPressure) THEN 00141 !For isothermal atmosphere, rho = rho_amb * exp(-phi/T_amb) 00142 !and pressure goes as P_amb = rho_amb * T_amb * exp(-phi/T_amb) 00143 !So we need to solve for rho_amb so that P_amb = P_disk at r=r_disk 00144 ! rho_amb = P_disk/T_amb/exp(-phi/T_amb) 00145 phi=GravityPotential(PointGravity%mass, (/disk_radius/lScale,0d0,0d0/), PointGravity%soft_length,PointGravity%soft_function) 00146 ambient%density=(Disk%pressure/(ambient_temperature/TempScale))/exp(-phi/(ambient_temperature/TempScale)) 00147 ! More Ambient stuff to get pressure matching at disk edge 00148 CALL CreateProfile(Ambient%profile, nPoints, (/Mass_Field, P_Field/), RADIAL) 00149 DO i=1,nPoints 00150 r=REAL(i-1)/REAL(nPoints)*r_max 00151 !Here we are not just using softening to get ambient density near star, but also fixing the ambient density within the sphere surrounding the disk to be constant to avoid Nan's in the ambient solution 00152 phi=GravityPotential(PointGravity%mass, (/max(r, disk_radius/lScale),0d0,0d0/), PointGravity%soft_length,PointGravity%soft_function) 00153 rho=ambient%density*exp(-phi/(ambient_temperature/TempScale)) 00154 Ambient%profile%data(i,1:2) = (/r, rho/) 00155 END DO 00156 00157 CALL Profile_PointGravityHSE(Ambient%profile, PointGravity, ambient%pressure) 00158 END IF 00159 00160 CALL UpdateAmbient(Ambient) 00161 00162 00163 ! Stream stuff 00164 stream_density=stream_density/rScale 00165 stream_temperature=stream_temperature/TempScale 00166 stream_velocity=stream_velocity*1e5/velscale 00167 stream_xpos=stream_xpos/lScale 00168 stream_halfwidth=stream_halfwidth/lScale 00169 stream_halfheight=stream_halfheight/lScale 00170 CALL AddTracer(stream_tracer, 'Stream Tracer') 00171 00172 ! Cooling Stuff 00173 IF (lCooling) THEN 00174 IF (.NOT. lRestart) THEN 00175 CALL CreateCoolingObject(coolingobj) 00176 coolingobj%iCooling = DMCool 00177 coolingobj%floortemp = MinTemp 00178 coolingobj%mintemp = MinTemp 00179 END IF 00180 END IF 00181 00182 ! Any diagnostics in chombo file? 00183 CALL AddDiagnosticVar(CoolingStrength_Field, 'Cooling Strength') 00184 00185 ! Refinement Stuff 00186 CALL ClearAllRefinements() 00187 00188 ! Refine Disk 00189 CALL CreateRefinement(DiskRefinement) 00190 CALL CreateShape(DiskRefinement%Shape, CYLINDER, 1.2*(/Disk%radius, Disk%radius, Disk%height/)) 00191 DiskRefinement%MaxLevel=disk_maxlevel 00192 00193 ! Derefine Disk center 00194 CALL CreateRefinement(DiskAntiRefinement) 00195 DiskAntiRefinement%Tolerance=DEREFINE_INSIDE 00196 DiskAntiRefinement%MaxLevel=disk_maxlevel 00197 CALL CreateShape(DiskAntiRefinement%Shape, CYLINDER, .5*(/Disk%radius, Disk%radius, Disk%height/)) 00198 00199 ! Hotspot refinement 00200 IF (lRefineHotSpot) THEN 00201 CALL CreateRefinement(HotSpotRefinement) 00202 CALL CreateShape(HotSpotRefinement%Shape, RECTANGULAR_PRISM, hotspot_width*(/half,half,half/)) 00203 HotSpotRefinement%Shape%position=hotspot_pos 00204 END IF 00205 00206 ! Stream refinement of injection region 00207 CALL CreateRefinement(StreamInjectionRefinement) 00208 CALL CreateShape(StreamInjectionRefinement%Shape, RECTANGULAR_PRISM, 1.8*stream_radius*(/stream_halfwidth, stream_halfwidth, stream_halfheight/)) !make it 20% larger than stream 00209 StreamInjectionRefinement%Shape%position=(/stream_xpos, GxBounds(2,1), 0d0/) 00210 00211 ! Stream refinement based on gradients in stream_density 00212 CALL CreateRefinement(StreamRefinement) 00213 StreamRefinement%Field=stream_tracer 00214 StreamRefinement%scale=LINEARSCALE 00215 StreamRefinement%tolerance=1d0/stream_density 00216 00217 ! Processing Stuff 00218 CALL CreateProjection(Projection) 00219 Projection%Field%id=Mass_Field 00220 Projection%Field%component=GASCOMP 00221 Projection%Field%name='ColumnDensity' 00222 Projection%dim=3d0 00223 00224 CALL CreateProjection(Projection) 00225 Projection%Field%id=Mass_Field 00226 Projection%Field%component=GASCOMP 00227 Projection%Field%name='ColumnDensity' 00228 Projection%dim=1d0 00229 00230 END SUBROUTINE ProblemModuleInit 00231 00234 SUBROUTINE ProblemGridInit(Info) 00235 TYPE(InfoDef) :: Info 00236 END SUBROUTINE ProblemGridInit 00237 00240 SUBROUTINE ProblemBeforeStep(Info) 00241 TYPE(InfoDef) :: Info 00242 INTEGER :: i,j,k, mbc(3) 00243 REAL(KIND=xPREC) :: pos(3), dx(3), rho, r, x, z 00244 00245 !determine number of ghost zones for each dimension 00246 mbc=levels(Info%level)%gmbc(levels(Info%level)%step)*merge((/1,1,1/),(/0,0,0/),nDim>=(/1,2,3/)) 00247 dx=levels(Info%level)%dx*merge((/1,1,1/),(/0,0,0/),nDim>=(/1,2,3/)) 00248 00249 ! Initialize stream from y=0 boundary 00250 DO j=1-mbc(2), Info%mX(2)+mbc(2) 00251 pos(2)=Info%xBounds(2,1)+(REAL(j, KIND=qPREC)-half)*dx(2) 00252 IF (pos(2) < GxBounds(2,1)) THEN 00253 DO i=1-mbc(1), Info%mX(1)+mbc(1) 00254 pos(1)=Info%xBounds(1,1)+(REAL(i, KIND=qPREC)-half)*dx(1) 00255 x = (pos(1)-stream_xpos)/stream_halfwidth 00256 DO k=1-mbc(3), Info%mX(3)+mbc(3) 00257 pos(3)=Info%xBounds(3,1)+(REAL(k, KIND=qPREC)-half)*dx(3) 00258 z = pos(3)/stream_halfheight 00259 r=sqrt(x**2+z**2) 00260 rho=stream_density*exp(-r**2) 00261 IF (rho >= ambient%density) THEN 00262 IF (stream_tracer /= 0) Info%q(i,j,k,stream_tracer)=rho 00263 Info%q(i,j,k,irho) = rho 00264 Info%q(i,j,k,ivx)=0d0 00265 Info%q(i,j,k,ivy)=rho*stream_velocity 00266 IF (ivz /= 0) Info%q(i,j,k,ivz)=0d0 00267 IF (iE /= 0) Info%q(i,j,k,iE)=gamma7*rho*stream_temperature+half*rho*stream_velocity**2 00268 END IF 00269 END DO 00270 END DO 00271 END IF 00272 END DO 00273 END SUBROUTINE ProblemBeforeStep 00274 00277 SUBROUTINE ProblemAfterStep(Info) 00278 TYPE(InfoDef) :: Info 00279 END SUBROUTINE ProblemAfterStep 00280 00283 SUBROUTINE ProblemSetErrFlag(Info) 00284 TYPE(InfoDef) :: Info 00285 END SUBROUTINE ProblemSetErrFlag 00286 00287 SUBROUTINE ProblemBeforeGlobalStep(n) 00288 INTEGER :: n 00289 END SUBROUTINE ProblemBeforeGlobalStep 00290 00291 END MODULE Problem 00292