Scrambler  1
BasicDisk/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 BasicDisk 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 Disks
00037   USE Winds
00038   USE DataDeclarations
00039   USE Ambients
00040   USE PointGravitySrc
00041   USE ParticleDeclarations !n
00042   USE Refinements
00043   IMPLICIT NONE
00044   SAVE
00045 
00046   PUBLIC ProblemModuleInit, ProblemGridInit, &
00047        ProblemBeforeStep, ProblemAfterStep, ProblemSetErrFlag, ProblemBeforeGlobalStep
00048   TYPE(DiskDef), POINTER :: mydisk
00049   TYPE(WindDef), POINTER :: wind
00050   INTEGER :: nWinds
00051   TYPE(AmbientDef), POINTER :: Ambient
00052   REAL(KIND=qPREC) :: iAccretion=KRUMHOLZ_ACCRETION !n
00053   REAL(KIND=qPREC) :: Height=1   !Disk peak density
00054   REAL(KIND=qPREC) :: soft_radius=1d0, radius=10d0
00055   REAL(KIND=qPREC), DIMENSION(3) :: xloc=(/0,0,0/)  !Disk location      
00056   TYPE(RefinementDef), POINTER :: Refinement
00057 CONTAINS
00058 
00060    SUBROUTINE ProblemModuleInit()
00061       INTEGER :: i,edge,j, soft_function
00062       TYPE(ParticleDef), POINTER :: Particle !n
00063       REAL(KIND=qPREC) :: thickness=0 !thickness of disk smoothing region
00064       REAL(KIND=qPREC) :: ddensity=10d0, adensity=1d0   !Disk peak density
00065       REAL(KIND=qPREC) :: dtemp=10d0, atemp=10d0
00066       REAL(KIND=qpREC) :: mass=1d0
00067       REAL(KIND=qPREC), DIMENSION(3) :: velocity= (/0,0,0/)  !Disk velocity (in direction of disk axis)
00068       REAL(KIND=qPREC) :: theta=0     !Angle between X-axis and disk axis (towards y-axis)
00069       REAL(KIND=qPREC) :: phi=0       !Angle around X-axis to disk axis
00070       LOGICAL :: lSinkParticle=.true.
00071       !     REAL(KIND=qPREC) :: mu=0        !Defines ratio of maximum magnetic pressure to ambient pressure
00072       !     REAL(KIND=qPREC) :: eta=0       !Parameter that determines ratio of maximum poloidal pressure to maximum toroidal pressure
00073       !     REAL(KIND=qPREC) :: B_theta=0   !Angle from disk axis to define disk orientation
00074       !     REAL(KIND=qPREC) :: B_phi=0     !Rotation around disk axis to define disk orientation(velocity)
00075 !      REAL(KIND=qPREC) :: B_tor=0     !Maximum Bfield for toroidal configuration
00076 !      REAL(KIND=qPREC) :: B_pol=0     !Maximum Bfield for poloidal configuration
00077 !      REAL(KIND=qPREC) :: Omega=0     !Solid body angular rotation
00078 !      REAL(KIND=qPREC) :: m2A         !Azimuthal density perturbation
00079 !      INTEGER :: iTracer                                !Disk Tracer
00080       
00081 !      INTEGER :: nwaves=0, nMHDwaves=0
00082 !      REAL(KIND=qPREC), DIMENSION(3) :: wavevector
00083 !      REAL(KIND=qPREC) :: phase, amplitude, amplitudes(3)
00084 
00085 !      REAL(KIND=qPREC) :: alpha        !< True love alpha parameter (ratio of thermal to gravitational energy). \f$ \alpha = \frac{5}{2} \left ( \frac{3}{4 \pi \rho_o M^2} \right )^(1/3) \frac{c_s^2}{G} \f$
00086 !      REAL(KIND=qPREC) :: beta_rot   !< True love beta rotational parameter (ratio of rotational to gravitational energy). \f$ \beta_{\Omega} = \frac{1}{4 \pi} \frac{\Omega^2}{G \rho_o} \f$
00087 
00088       TYPE(PointGravityDef), POINTER :: PointGravityObj
00089 
00090       NAMELIST /ProblemData/ ddensity, adensity, velocity, xloc, radius, thickness, dtemp, atemp, theta, phi, mass, soft_radius, soft_function, height, lSinkParticle
00091 
00092       OPEN(UNIT=PROBLEM_DATA_HANDLE, FILE='problem.data', STATUS="OLD")
00093       READ(PROBLEM_DATA_HANDLE,NML=ProblemData)
00094       CALL CreateAmbient(Ambient)
00095       Ambient%density=adensity/rScale !cu
00096       Ambient%pressure= adensity*atemp/rScale/TempScale !cu
00097 
00098       NULLIFY(PointGravityObj)
00099       IF (lSinkParticle) THEN
00100          NULLIFY(Particle) !n
00101          CALL CreateParticle(Particle) !n
00102          Particle%q(1)=mass*mSolar/mScale !n
00103          Particle%xloc=xloc !n
00104          Particle%iAccrete=iAccretion !n
00105          CALL CreatePointGravityObject(Particle%PointGravityObj)
00106          Particle%lFixed=.true. !n
00107          CALL AddSinkParticle(Particle)
00108          PointGravityObj=>Particle%PointGravityObj
00109       ELSE
00110          CALL CreatePointGravityObject(PointGravityObj)
00111       END IF
00112    
00113          !        CALL CreatePointGravityObject(PointGravityObj)    
00114       PointGravityObj%mass=mass*mSolar/mScale
00115       PointGravityObj%x0=xloc
00116       PointGravityObj%v0=velocity
00117       PointGravityObj%t0=levels(0)%tnow
00118       PointGravityObj%soft_length=soft_radius*radius/lScale
00119       PointGravityObj%soft_function=soft_function
00120 
00121       CALL CreateDisk(mydisk)
00122 
00123       mydisk%soft_length=soft_radius*radius/lScale
00124       mydisk%soft_function=soft_function
00125       mydisk%pressure=ddensity*dtemp/rScale/TempScale !cu
00126       mydisk%density=ddensity/rScale
00127       !mydisk%Type=HYDROSTATIC
00128       ! Calculate solid body rotational velocity by beta_rot
00129 !      mydisk%omega = sqrt(4d0*pi*ScaleGrav*density*beta_rot)
00130       mydisk%velocity=velocity
00131       mydisk%theta=theta
00132       mydisk%phi=phi
00133 !      mydisk%B_tor=B_tor
00134 !      mydisk%B_pol=B_pol
00135       mydisk%position=xloc
00136       mydisk%thickness=4d0*sink_dx
00137       mydisk%radius=radius/lScale
00138       mydisk%height=height/lScale
00139       mydisk%central_mass=mass*mSolar/mScale
00140       !mydisk%Background_density=ambient%density
00141       !mydisk%SubSample=0
00142 !      mydisk%m2A=m2A
00143       CALL UpdateDisk(mydisk)
00144       
00145       CLOSE(PROBLEM_DATA_HANDLE)
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%Type=OUTFLOW_ONLY
00154                Wind%dir=i
00155                Wind%edge=edge
00156             END IF
00157          END DO
00158       END DO
00159 
00160 
00161       CALL CreateRefinement(Refinement)
00162       CALL CreateShape(Refinement%Shape)
00163       CALL SetShapeType(Refinement%Shape, CYLINDER, (/myDisk%radius,myDisk%radius, myDisk%Height/))
00164       Refinement%BufferCells=4
00165 
00166    END SUBROUTINE ProblemModuleInit
00167 
00170   SUBROUTINE ProblemGridInit(Info)
00171     TYPE(InfoDef) :: Info
00172   END SUBROUTINE ProblemGridInit
00173 
00176   SUBROUTINE ProblemBeforeStep(Info)
00177     TYPE(InfoDef) :: Info
00178     INTEGER :: i
00179 !    DO i=1,nWinds
00180 !       CALL BeforeStepWind(Info,Wind(i))
00181 !    END DO
00182   END SUBROUTINE ProblemBeforeStep
00183 
00186   SUBROUTINE ProblemAfterStep(Info)
00187     TYPE(InfoDef) :: Info
00188   END SUBROUTINE ProblemAfterStep
00189 
00192   SUBROUTINE ProblemSetErrFlag(Info)
00193     TYPE(InfoDef) :: Info
00194       INTEGER :: i,j,k,mx,my,mz,rmbc,zrmbc,level
00195       REAL(KIND=xprec) :: x,y,z,xl,yl,zl,dx,dy,dz
00196 
00197 
00198       level=Info%level
00199       rmbc=levels(level)%gmbc(levels(level)%step)
00200       zrmbc=rmbc
00201                 mx = Info%mX(1)
00202                 my = Info%mX(2)
00203                 mz = Info%mX(3)
00204       dx=levels(level)%dX
00205       dy=dx;dz=dx
00206       xl=Info%xBounds(1,1)
00207       yl=Info%xBounds(2,1)
00208       zl=Info%xBounds(3,1)
00209 
00210       DO k=1, mz ; z = (zl+(REAL(k,xPrec)-half)*dz)
00211                    if (abs(z-xloc(3)).gt.height/lScale+2d0*dz) cycle 
00212               DO i=1, mx ; x = (xl+(REAL(i,xPrec)-half)*dx)
00213             DO j=1, my ; y = (yl+(REAL(j,xPrec)-half)*dy)
00214                          if (sqrt((x-xloc(1))**2+(y-xloc(2))**2).le.radius/lScale+2d0*dx) Info%ErrFlag(i,j,k)=1
00215       END DO; END DO; END DO
00216   END SUBROUTINE ProblemSetErrFlag
00217 
00218   SUBROUTINE ProblemBeforeGlobalStep(n)
00219      INTEGER :: n
00220   END SUBROUTINE ProblemBeforeGlobalStep
00221 
00222 END MODULE Problem
00223 
 All Classes Files Functions Variables