Scrambler  1
ThermalInstability/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 ThermalInstability 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 Clumps
00037   USE Winds
00038   USE DataDeclarations
00039   USE CoolingSrc
00040   USE Ambients
00041   IMPLICIT NONE
00042   SAVE
00043   PUBLIC ProblemModuleInit, ProblemGridInit, &
00044        ProblemBeforeStep, ProblemAfterStep, ProblemSetErrFlag, ProblemBeforeGlobalStep
00045   TYPE(ClumpDef), POINTER :: myclump
00046   TYPE(WindDef), POINTER  :: Wind
00047   INTEGER :: nWinds
00048   TYPE(CoolingDef),POINTER :: coolingobj
00049 CONTAINS
00050 
00052    SUBROUTINE ProblemModuleInit()
00053       INTEGER :: i,edge,j
00054       REAL(KIND=qPREC) :: radius=1    !clump radius
00055       REAL(KIND=qPREC) :: thickness=.1 !thickness of clump smoothing region
00056       REAL(KIND=qPREC) :: density=10   !Clump peak density
00057       REAL(KIND=qPREC) :: temp=.1      !Clump temperature
00058       REAL(KIND=qPREC), DIMENSION(3) :: velocity= (/0,0,0/)  !Clump velocity (in direction of clump axis)
00059       REAL(KIND=qPREC) :: theta=0     !Angle between X-axis and clump axis (towards y-axis)
00060       REAL(KIND=qPREC) :: phi=0       !Angle around X-axis to clump axis
00061       !     REAL(KIND=qPREC) :: mu=0        !Defines ratio of maximum magnetic pressure to ambient pressure
00062       !     REAL(KIND=qPREC) :: eta=0       !Parameter that determines ratio of maximum poloidal pressure to maximum toroidal pressure
00063       !     REAL(KIND=qPREC) :: B_theta=0   !Angle from clump axis to define clump orientation
00064       !     REAL(KIND=qPREC) :: B_phi=0     !Rotation around clump axis to define clump orientation(velocity)
00065       REAL(KIND=qPREC) :: B_tor=0     !Maximum Bfield for toroidal configuration
00066       REAL(KIND=qPREC) :: B_pol=0     !Maximum Bfield for poloidal configuration
00067       REAL(KIND=qPREC) :: Omega=0     !Solid body angular rotation
00068       REAL(KIND=qPREC) :: m2A         !Azimuthal density perturbation
00069       REAL(KIND=qPREC), DIMENSION(3) :: xloc=(/0,0,0/)  !Clump location
00070       INTEGER :: iTracer                                !Clump Tracer
00071 
00072       INTEGER :: nwaves=0, nMHDwaves=0
00073       REAL(KIND=qPREC), DIMENSION(3) :: wavevector
00074       REAL(KIND=qPREC) :: phase, amplitude, amplitudes(3)
00075 
00076       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$
00077       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$
00078       LOGICAL :: lCooling=.false.
00079       TYPE(AmbientDef), POINTER :: Ambient
00080       REAL(KIND=qPREC) :: rhoOut, pOut, vxOut, vyOut, vzOut, BxOut, ByOut, BzOut
00081       NAMELIST /AmbientData/ rhoOut, pOut, vxOut, vyOut, vzOut, BxOut, ByOut, BzOut
00082       NAMELIST /ProblemData/ density,velocity, xloc, radius, thickness, temp,theta, phi, B_tor, B_pol, beta_rot, alpha, m2A, nwaves, nMHDwaves, lCooling
00083       NAMELIST /WaveData/ wavevector, amplitude, phase
00084       NAMELIST /MHDWaveData/ wavevector, amplitudes, phase
00085 
00086       OPEN(UNIT=PROBLEM_DATA_HANDLE, FILE='problem.data', STATUS="OLD")
00087       READ(PROBLEM_DATA_HANDLE,NML=ProblemData)
00088       READ(PROBLEM_DATA_HANDLE,NML=AmbientData)
00089 
00090       CALL CreateAmbient(Ambient)
00091       Ambient%density=rhoOut
00092       Ambient%pressure=pOut
00093       Ambient%B(:)=(/BxOut, ByOut, BzOut/)
00094       Ambient%velocity(:)=(/vxOut, vyOut, vzOut/)
00095 
00096       IF (lCooling) THEN
00097          IF (.NOT. lRestart) THEN
00098             ! see sources/cooling.f90::CreateCoolingObject for
00099             ! default values of a cooling source term
00100             CALL CreateCoolingObject(coolingobj)
00101          ELSE
00102             coolingobj => firstcoolingobj
00103          END IF
00104          coolingobj%iCooling=IICOOL
00105          coolingobj%floortemp=1d0
00106          coolingobj%mintemp=0.001
00107       END IF
00108 
00109 
00110       CALL CreateClump(myclump)
00111       ! Modify density or temperature based on alpha parameter
00112       IF (alpha == 0 .OR. .NOT. lSelfGravity) THEN
00113          !         myclump%density=density
00114          IF (lSelfGravity) write(*,*) 'alpha = ', 5d0/2d0*gamma*Temp/ScaleGrav/density/(4d0/3d0*pi*radius**2)
00115          IF (temp == 0) THEN
00116             CALL InitIICool(coolingobj)
00117             temp=GetIICoolEqTemp(density)
00118             Ambient%Pressure=Ambient%density*temp/TempScale
00119             write(*,*) 'setting Ambient%Pressure=', Ambient%Pressure
00120          END IF
00121       ELSE
00122          IF (iEOS == EOS_ISOTHERMAL) THEN
00123             density=5d0/2d0*Iso_Speed2/ScaleGrav/alpha/(4d0/3d0*pi*radius**2)
00124             write(*,*) "adjusting clump density to", density
00125          ELSE
00126             temp=2d0/5d0*density*ScaleGrav*alpha*(4d0/3d0*pi*radius**2)/gamma
00127             write(*,*) "adjusting clump temperature to", temp
00128          END IF
00129       END IF
00130 
00131 
00132       myclump%temperature=0
00133       myclump%pressure=density*temp
00134       myclump%density=density
00135       ! Calculate solid body rotational velocity by beta_rot
00136       myclump%omega = sqrt(4d0*pi*ScaleGrav*density*beta_rot)
00137       myclump%velocity=velocity
00138       myclump%theta=theta
00139       myclump%phi=phi
00140       myclump%B_toroidal=B_tor
00141       myclump%B_poloidal=B_pol
00142       myclump%position=xloc
00143       myclump%thickness=thickness
00144       myclump%radius=radius
00145       myclump%m2A=m2A
00146       CALL UpdateClump(myclump)
00147 
00148       IF (nwaves > 0) THEN
00149          ALLOCATE(myclump%DensityPerturbation)
00150          myclump%DensityPerturbation%type=COSINESERIES
00151          myclump%DensityPerturbation%geometry=SPHERICAL
00152          CALL InitPerturbationWaves(Myclump%DensityPerturbation, nwaves)
00153          DO i=1,nWaves
00154             READ(PROBLEM_DATA_HANDLE, NML=WaveData)
00155             CALL AddPerturbationWave(Myclump%DensityPerturbation, wavevector, phase, amplitude)
00156          END DO
00157       END IF
00158 
00159       IF (nMHDwaves > 0) THEN
00160          ALLOCATE(myclump%MagneticPerturbation(nEmf))        
00161          DO i=1,nEMF
00162             myclump%MagneticPerturbation%type=COSINESERIES
00163             CALL InitPerturbationWaves(Myclump%MagneticPerturbation(i), nMHDwaves)
00164          END DO
00165          DO i=1,nMHDWaves
00166             READ(PROBLEM_DATA_HANDLE, NML=MHDWaveData)
00167             DO j=1,nEMF
00168                CALL AddPerturbationWave(Myclump%MagneticPerturbation(j), wavevector, phase, amplitudes(j))
00169             END DO
00170          END DO
00171       END IF
00172 
00173 
00174       CLOSE(PROBLEM_DATA_HANDLE)
00175 
00176       nWinds=0
00177       DO i=1,nDim
00178          DO edge=1,2
00179             IF (Gmthbc(i,edge) == 1) THEN 
00180                nWinds=nWinds+1
00181                CALL CreateWind(Wind)
00182                Wind%dir=i
00183                Wind%edge=edge
00184                Wind%density=Ambient%density
00185                Wind%Temperature=temp
00186             END IF
00187          END DO
00188       END DO
00189 
00190    END SUBROUTINE ProblemModuleInit
00191 
00194   SUBROUTINE ProblemGridInit(Info)
00195     TYPE(InfoDef) :: Info
00196   END SUBROUTINE ProblemGridInit
00197 
00200   SUBROUTINE ProblemBeforeStep(Info)
00201     TYPE(InfoDef) :: Info
00202     INTEGER :: i
00203 !    DO i=1,nWinds
00204 !       CALL BeforeStepWind(Info,Wind(i))
00205 !    END DO
00206   END SUBROUTINE ProblemBeforeStep
00207 
00210   SUBROUTINE ProblemAfterStep(Info)
00211     TYPE(InfoDef) :: Info
00212   END SUBROUTINE ProblemAfterStep
00213 
00216   SUBROUTINE ProblemSetErrFlag(Info)
00217     TYPE(InfoDef) :: Info
00218   END SUBROUTINE ProblemSetErrFlag
00219 
00220   SUBROUTINE ProblemBeforeGlobalStep(n)
00221      INTEGER :: n
00222   END SUBROUTINE ProblemBeforeGlobalStep
00223 
00224 END MODULE Problem
00225 
 All Classes Files Functions Variables