Scrambler  1
disks.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 !    disks.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 
00029 
00032 
00036 
00039 MODULE Disks
00040    USE GlobalDeclarations
00041    USE DataDeclarations
00042    USE PhysicsDeclarations
00043    USE CommonFunctions
00044    USE Shapes
00045    USE Perturbation
00046    USE DataInfoOps
00047    USE EOS
00048    USE BE_MODULE
00049    USE ObjectDeclarations
00050    IMPLICIT NONE
00051 
00052    INTEGER, PARAMETER :: UNIFORM = 0, HYDROSTATIC = 1
00053 
00055    TYPE DiskDef
00056       REAL(KIND=qPREC), DIMENSION(3) :: velocity= (/0,0,0/)  !Disk velocity 
00057       REAL(KIND=qPREC), DIMENSION(3) :: position = (/0,0,0/)  !Disk velocity 
00058       REAL(KIND=qPREC) :: density=10   !Disk peak density
00059       REAL(KIND=qPREC) :: pressure=1   !External Disk Pressure
00060       REAL(KIND=qPREC) :: radius=1
00061       REAL(KIND=qPREC) :: height=.25
00062       REAL(KIND=qPREC) :: thickness=.1 !thickness of disk smoothing region
00063       REAL(KIND=qPREC) :: phi=0
00064       REAL(KIND=qPREC) :: theta=0
00065       REAL(KIND=qPREC) :: central_mass  !Central Mass
00066 !      REAL(KIND=qPREC) :: temp=.1      !Disk temperature
00067       REAL(KIND=qPREC) :: t0=0         !Last time velocity and position were updated.
00068       REAL(KIND=qPREC) :: background_density=1   !Last time velocity and position were updated.
00069 !      REAL(KIND=qPREC) :: alpha=.5d0
00070 !      REAL(KIND=qPREC) :: B_tor=0     !Maximum Bfield for toroidal configuration
00071 !      REAL(KIND=qPREC) :: B_pol=0     !Maximum Bfield for poloidal configuration
00072 !      REAL(KIND=qPREC) :: Omega=0     !Angular velocity for solid body rotation
00073       INTEGER :: soft_function         !Gravitational softening function
00074       REAL(KIND=qPREC) :: soft_length  !Gravitational softening radius
00075       INTEGER :: SubSample=1           !Subsample factor
00076       INTEGER :: iTracer=0             !Disk Tracer
00077       INTEGER :: type = UNIFORM
00078       LOGICAL :: PersistInternal=.false.
00079       LOGICAL :: PersistInBoundaries(3,2)=.false.  !keep disk in each boundary
00080 !     REAL(KIND=qPREC) :: m2A=0 !Amplitude of an m=2 density perturbation [  rho = rho_0*(1+m2A*cos(2*phi))   ]
00081       INTEGER :: id
00082       TYPE(ShapeDef), POINTER :: Shape
00083 !    TYPE(PerturbationDef), POINTER :: DensityPerturbation => Null()
00084 !    TYPE(PerturbationDef), DIMENSION(:), POINTER :: MagneticPerturbation
00085       INTEGER :: ObjId
00086    END TYPE DiskDef
00087 
00088   !new declaration
00089   TYPE pDiskDef
00090     TYPE(DiskDef), POINTER :: ptr
00091   END TYPE pDiskDef
00092   TYPE(pDiskDef) :: pDisk
00093   !
00094 
00095    SAVE
00096 CONTAINS
00097 
00098 
00101    SUBROUTINE CreateDisk(Disk, density, pressure, central_mass)
00102       TYPE(DiskDef), POINTER :: Disk
00103       REAL(KIND=qPREC), OPTIONAL :: density, pressure, central_mass
00104       IF (.NOT. ASSOCIATED(Disk)) ALLOCATE(Disk)  
00105       ALLOCATE(disk%Shape)
00106       IF(Present(density)) Disk%density = density
00107       IF(Present(pressure)) Disk%pressure = pressure
00108       IF(Present(central_mass)) Disk%central_mass = central_mass 
00109       CALL UpdateDisk(Disk)
00110       CALL AddDiskToList(Disk)
00111    END SUBROUTINE CreateDisk
00112 
00115    SUBROUTINE UpdateDisk(Disk)
00116       TYPE(DiskDef), POINTER :: Disk
00117       IF (Disk%Type == UNIFORM) THEN
00118          CALL SetShapeType(Disk%Shape, CYLINDER, (/Disk%radius,Disk%radius, Disk%height/))
00119       ELSE
00120          CALL SetShapeType(Disk%Shape, CYLINDER, (/Disk%radius*1e2,Disk%radius*1e2, Disk%height*1e2/))
00121       END IF
00122       CALL SetShapeOrientation(Disk%Shape, 0d0, Disk%theta, Disk%phi)
00123       Disk%Shape%Position=Disk%position
00124       CALL SetShapeBounds(Disk%Shape)         
00125       
00126    END SUBROUTINE UpdateDisk
00127 
00128    SUBROUTINE AddDiskToList(Disk)
00129       TYPE(DiskDef), POINTER :: Disk
00130       TYPE(ObjectDef), POINTER :: Object
00131       Disk%ObjId = ObjectListAdd(Object,DISKOBJ)
00132       pDisk%ptr => Disk
00133       len = size(transfer(pDisk, dummy_char))
00134       ALLOCATE(Object%storage(len))
00135       Object%storage = transfer(pDisk, Object%storage)
00136    END SUBROUTINE AddDiskToList
00137 
00138    SUBROUTINE DestroyDisk(Disk,id)
00139       TYPE(DiskDef),POINTER :: Disk
00140       TYPE(ObjectDef),POINTER :: Object
00141       INTEGER,OPTIONAL :: id
00142 
00143       IF(PRESENT(id)) THEN
00144         Object => ObjectListFind(id)
00145         IF (ASSOCIATED(Object) .AND. Object%type == DISKOBJ) THEN
00146           write(*,*) "found obj"
00147           pDisk = transfer(Object%storage,pDisk)
00148           IF(ASSOCIATED(pDisk%ptr%Shape)) THEN
00149             DEALLOCATE(pDisk%ptr%Shape)
00150           ENDIF
00151           DEALLOCATE(pDisk%ptr)
00152           NULLIFY(pDisk%ptr)
00153           CALL ObjectListRemove(id)
00154         ENDIF
00155       ELSE
00156         CALL ObjectListRemove(Disk%ObjId)
00157         IF(ASSOCIATED(Disk%Shape)) THEN
00158             DEALLOCATE(Disk%Shape)
00159           ENDIF
00160         DEALLOCATE(Disk)
00161         NULLIFY(Disk)
00162       ENDIF
00163    END SUBROUTINE DestroyDisk
00164 
00165    SUBROUTINE DiskGridInit(Info, Disk)
00166       TYPE(InfoDef) :: Info
00167       TYPE(DiskDef), POINTER :: Disk
00168       INTEGER, POINTER, DIMENSION(:,:,:) :: mSs
00169       REAL(KIND=qPREC), POINTER, DIMENSION(:,:) :: offsets
00170       INTEGER :: nOverlaps
00171       CALL CalcPhysicalOverlaps(Info, Disk%Shape%xBounds, mSs, nOverlaps, offsets, IEVERYWHERE, lHydroPeriodic)
00172       IF (nOverlaps > 0) THEN
00173          CALL PlaceDisk(Info, Disk, nOverlaps, mSs, offsets)
00174          DEALLOCATE(mSs, offsets)
00175       END IF
00176     END SUBROUTINE DiskGridInit
00177 
00178    
00179    SUBROUTINE DiskBeforeStep(Info, Disk)
00180       TYPE(InfoDef) :: Info
00181       TYPE(DiskDef), POINTER :: Disk
00182       INTEGER, POINTER, DIMENSION(:,:,:) :: mSs
00183       REAL(KIND=qPREC), POINTER, DIMENSION(:,:) :: offsets
00184       INTEGER :: nOverlaps
00185       INTEGER :: i,j
00186       IF (Disk%PersistInternal) THEN
00187          CALL DiskGridInit(Info, Disk)
00188       ELSE 
00189          DO i=1,nDim
00190             DO j=1,2
00191                IF (Disk%PersistInBoundaries(i,j)) THEN
00192                   CALL CalcPhysicalOverlaps(Info, Disk%Shape%xBounds, mSs, nOverlaps, offsets, IBOUNDARY(i,j), lHydroPeriodic)
00193                   IF (nOverlaps > 0) THEN
00194                      CALL PlaceDisk(Info, Disk, nOverlaps, mSs, offsets)
00195                      DEALLOCATE(mSs, offsets)
00196                   END IF
00197                END IF
00198             END DO
00199          END DO
00200       END IF
00201    END SUBROUTINE DiskBeforeStep
00202 
00203    SUBROUTINE DiskSetErrFlag(Info, Disk)
00204       TYPE(InfoDef) :: Info
00205       Type(DiskDef), POINTER :: Disk
00206    END SUBROUTINE DiskSetErrFlag
00207 
00208    SUBROUTINE DiskBeforeGlobalStep(n)
00209       INTEGER :: n
00210    END SUBROUTINE DiskBeforeGlobalStep
00211 
00215    SUBROUTINE PlaceDisk(Info,Disk,nOverlaps,mSs,offsets)
00216       TYPE(InfoDef) :: Info
00217       Type(DiskDef) :: Disk
00218       INTEGER :: i,j,k,n,m,ii,jj,kk, location
00219       INTEGER, DIMENSION(3,2) :: mS    
00220       INTEGER, POINTER, DIMENSION(:,:,:) :: mSs
00221       REAL(KIND=qPREC), DIMENSION(3) :: offset
00222       REAL(KIND=qPREC), POINTER, DIMENSION(:,:) :: offsets
00223       REAL(KIND=qPREC), DIMENSION(3) :: xpos, pos, rpos
00224       REAL(KIND=qPREC) :: sample_fact(3), q_fact, dx(3)
00225       INTEGER :: sample_res(3), nOverlaps
00226       REAL(KIND=qPREC), DIMENSION(3,2) :: tempbounds
00227       REAL(KIND=qPREC), DIMENSION(:), ALLOCATABLE :: q_source, emf_source, qacc
00228       REAL(KIND=qPREC), DIMENSION(:,:,:,:), ALLOCATABLE :: emf
00229 !      RETURN
00230       dx=levels(Info%level)%dX*(/(merge(1,0,nDim >= i), i=1,3)/)
00231       xpos=0          
00232       ALLOCATE(q_Source(NrHydroVars), qacc(NrHydroVars))
00233       sample_res=1
00234       sample_fact=0d0
00235       IF (Disk%SubSample > 0) sample_res(1:nDim)=max(Disk%SubSample,2**(MaxLevel-Info%level))
00236       sample_fact(1:nDim)=1d0/REAL(sample_res(1:nDim),8)
00237       q_fact=product(sample_fact(1:nDim))
00238       DO n=1,nOverlaps
00239          mS=mSs(n,:,:)
00240          offset=offsets(n,:)
00241          !Set up emf first - aux fields first
00242          CALL ConvertTotalToInternalEnergy(Info%q(mS(1,1):mS(1,2), mS(2,1):mS(2,2),mS(3,1):mS(3,2),:))            ! Now set up cell centered quantities (density and momentum)
00243             DO k=mS(3,1),mS(3,2)
00244                xpos(3)=Info%xBounds(3,1)+offset(3)+(k-1)*dx(3)
00245                DO j=mS(2,1),mS(2,2)
00246                   xpos(2)=Info%xBounds(2,1)+offset(2)+(j-1)*dx(2)
00247                   DO i=mS(1,1),mS(1,2)
00248                      xpos(1)=Info%xBounds(1,1)+offset(1)+(i-1)*dx(1)
00249                      qacc=0
00250                      DO kk=1,sample_res(3)
00251                         pos(3)=xpos(3)+(REAL(kk, 8)-half)*dx(3)*sample_fact(3)
00252                         DO jj=1,sample_res(2)
00253                            pos(2)=xpos(2)+(REAL(jj, 8)-half)*dx(2)*sample_fact(2)
00254                            DO ii=1,sample_res(1)
00255                               pos(1)=xpos(1)+(REAL(ii, 8)-half)*dx(1)*sample_fact(1)
00256                               IF (IsInShape(Disk%Shape, pos, rpos)) THEN
00257                                  CALL q_DiskSource(Disk,rpos,Info%q(i,j,k,1:NrHydroVars), q_source)
00258                                  qacc=qacc+q_source
00259 !                              ELSE
00260 !                                 write(*,*) pos
00261                               END IF
00262                            END DO
00263                         END DO
00264                      END DO
00265                      Info%q(i,j,k,1:NrHydroVars)=Info%q(i,j,k,1:NrHydroVars)+qacc*q_fact                 
00266                   END DO
00267                END DO
00268             END DO
00269             CALL ConvertInternalToTotalEnergy(Info%q(mS(1,1):mS(1,2), mS(2,1):mS(2,2),mS(3,1):mS(3,2),:))           
00270          END DO
00271          DEALLOCATE(q_Source, qacc)
00272    END SUBROUTINE PlaceDisk
00273 
00274    SUBROUTINE q_DiskSource(Disk, pos, q, s)
00275       TYPE(DiskDef) :: Disk
00276       REAL(KIND=qPREC), DIMENSION(:) :: pos, q, s
00277       REAL(KIND=qPREC) :: r, f, rotation_speed(3), perturb, r_s, omega, f_rad, f_grav(3), c_s, rho, r_disk
00278       !     r=SQRT(SUM((pos(1:nDim)-Disk%position(1:nDim))**2))
00279       r=sqrt(sum(pos**2))
00280       f=smooth_tanh(r/Disk%Shape%size_param(1),Disk%thickness)
00281       perturb=1d0
00282       s=0
00283       r_s=sqrt(sum(pos(1:2)**2))
00284 
00285       IF (Disk%Type == UNIFORM) THEN
00286          rho=Disk%density*perturb
00287          s(1)=(rho-q(1))*f !fade into disk        
00288          f_grav=GravityForce(Disk%central_mass, TransformCoordsFromShape(Disk%Shape,pos), Disk%soft_length, Disk%soft_function)
00289 
00290          !gravitational force of rotated position relative to particl
00291          f_grav=RotateVectorToShape(Disk%Shape, f_grav)
00292          IF (r > TINY(r)) THEN
00293             f_rad=sum(-f_grav(1:2)*pos(1:2))/r_s
00294             omega=sqrt(f_rad/r_s)
00295          ELSE
00296             omega=0
00297          END IF
00298       ELSEIF (Disk%Type == HYDROSTATIC) THEN
00299          IF (iEOS == EOS_ISOTHERMAL) THEN
00300             c_s=Iso_Speed
00301          ELSE
00302             c_s=sqrt(gamma*Disk%Pressure/Disk%density)
00303          END IF
00304          r_disk=Disk%radius
00305          IF (Disk%soft_function == PLUMMERSOFT) THEN
00306             r_s=sqrt(r_s**2+disk%soft_length**2)
00307             r=sqrt(r**2+disk%soft_length**2)
00308             r_disk=sqrt(r_disk**2+disk%soft_length**2)
00309          END IF
00310 
00311          IF (r_s > r_disk) THEN
00312             rho=Disk%density*perturb*exp(Disk%central_mass*ScaleGrav/c_s**2*(1d0/r-1d0/r_disk))
00313             omega=0
00314             IF (rho < Disk%Background_density) THEN
00315                rho=q(1)
00316                s(1)=0
00317             ELSE
00318                s(1)=(rho-q(1))*f
00319             END IF
00320          ELSE
00321             rho=Disk%density*perturb*exp(Disk%central_mass*ScaleGrav/c_s**2*(1d0/r-1d0/r_s))
00322             IF (rho < Disk%background_density) THEN
00323                rho=q(1)
00324                s(1)=0d0!(rho-q(1))*f
00325                omega=0d0!sqrt((Disk%alpha)*Disk%central_mass*ScaleGrav/r_s**2)
00326             ELSE
00327                s(1)=(rho-q(1))*f
00328                omega=sqrt(Disk%central_mass*ScaleGrav/r_s**3)
00329             END IF
00330 
00331         END IF
00332       END IF
00333       rotation_speed(1:3)=RotateVectorFromShape(Disk%Shape, Omega*(/-pos(2),pos(1),0d0/)) 
00334       s(ivx:ivx+nDim-1)=(rho*(Disk%velocity(1:nDim) + rotation_speed(1:nDim))-q(ivx:ivx+nDim-1))*f!momentum associated with disk mass + rotational energy
00335       IF (iCylindrical == WITHANGMOM) s(iAngMom)=(rho*(rotation_speed(3))-q(iAngMom))*f
00336       IF (iE .ne. 0) s(iE)=(gamma7*Disk%pressure-q(iE))*f !Internal energy associated with disk material
00337       IF (Disk%iTracer .ne. 0) s(Disk%iTracer) = Disk%density*perturb*f
00338    END SUBROUTINE q_DiskSource
00339 
00340 
00341 END MODULE Disks
00342 
00343 
 All Classes Files Functions Variables