Scrambler  1
ambients.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 !    ambients.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 MODULE Ambients
00033    USE GlobalDeclarations
00034    USE DataDeclarations
00035    USE PhysicsDeclarations
00036    USE EOS
00037    USE Perturbation
00038    USE VectorPerturbation
00039    USE ObjectDeclarations
00040    USE Profiles
00041    IMPLICIT NONE
00043    TYPE AmbientDef
00044       REAL(KIND=qPREC) :: density=1d0
00045       REAL(KIND=qPREC) :: velocity(3)=0d0
00046       REAL(KIND=qPREC) :: pressure=1d0
00047       REAL(KIND=qPREC) :: B(3) = 0d0
00048       REAL(KIND=qPREC), DIMENSION(3,2) :: xBounds
00049       TYPE(PerturbationDef), POINTER :: DensityPerturbation => null()
00050       TYPE(VectorPerturbationDef), POINTER :: VelocityPerturbation => null()
00051       TYPE(VectorPerturbationDef), POINTER :: MagneticPerturbation => null()
00052       TYPE(ProfileDef), POINTER :: profile => null()
00053       LOGICAL :: PersistInBoundaries(3,2) = .false.
00054       INTEGER :: iTracer=0
00055       INTEGER :: ObjId
00056       TYPE(PointGravityDef), POINTER :: PointGravity => null()
00057    END TYPE AmbientDef
00058 
00059    !new declaration
00060    TYPE pAmbientDef
00061       TYPE(AmbientDef), POINTER :: ptr
00062    END TYPE pAmbientDef
00063    TYPE(pAmbientDef) :: pAmbient
00064    !new declaration
00065 
00066 CONTAINS
00067 
00068 
00071    SUBROUTINE CreateAmbient(Ambient, density, pressure)
00072       TYPE(AmbientDef), POINTER :: Ambient
00073       REAL(KIND=qPREC), OPTIONAL :: density, pressure
00074       INTEGER :: i
00075       ALLOCATE(Ambient)
00076       IF (Present(density)) Ambient%density=density
00077       IF (Present(pressure)) Ambient%pressure=pressure
00078       Ambient%xBounds(:,1)=GxBounds(:,1)-(/(merge(1,0,ndim>=i),i=1,3)/)*levels(0)%dx*levels(0)%gmbc(1)
00079       Ambient%xBounds(:,2)=GxBounds(:,2)+(/(merge(1,0,ndim>=i),i=1,3)/)*levels(0)%dx*levels(0)%gmbc(1)
00080       CALL AddAmbientToList(Ambient)
00081    END SUBROUTINE CreateAmbient
00082 
00083 
00086    SUBROUTINE UpdateAmbient(Ambient)
00087      TYPE(AmbientDef), POINTER :: Ambient
00088      INTEGER :: i,j
00089      !Do a sanity check on magnetic fields and boundary conditions
00090      DO i=1,nDim
00091         IF (ANY(Gmthbc(i,:) == REFLECT_WALL) .AND. Ambient%B(i) /= 0) THEN
00092            print*, 'cannot have boundary conditions of type ', REFLECT_WALL, ' along dimension ', i, ' with B(',i,') /= 0'
00093            STOP
00094         END IF
00095         
00096         IF (ANY(Gmthbc(i,:) == REFLECT_BPARALLEL) .AND. ANY(Ambient%B(modulo(i+(/(j,j=0,nDim-2)/),nDim)+1) /= 0)) THEN
00097            print*, 'cannot have boundary conditions of type ', REFLECT_BPARALLEL,' along dimension ', i, &
00098                 ' with B(',modulo(i,3)+1,') or B(',modulo(i+1,3)+1,') /= 0'
00099            STOP
00100         END IF
00101      END DO
00102    END SUBROUTINE UpdateAmbient
00103 
00104 
00105    SUBROUTINE AddAmbientToList(Ambient)
00106      TYPE(AmbientDef), POINTER :: Ambient
00107      TYPE(ObjectDef), POINTER :: Object
00108      Ambient%ObjId = ObjectListAdd(Object,AMBIENTOBJ)
00109      pAmbient%ptr => Ambient
00110      len = size(transfer(pAmbient, dummy_char))
00111      ALLOCATE(Object%storage(len))
00112      Object%storage = transfer(pAmbient, Object%storage)
00113    END SUBROUTiNE AddAmbientToList
00114 
00115    SUBROUTINE DestroyAmbient(Ambient)
00116       TYPE(AmbientDef), POINTER :: Ambient
00117       TYPE(ObjectDef), POINTER :: Object
00118       INTEGER :: id
00119       id=Ambient%ObjId
00120       Object => ObjectListFind(id)
00121       IF (ASSOCIATED(Object) .AND. Object%type == AMBIENTOBJ) THEN
00122          pAmbient = transfer(Object%storage,pAmbient)
00123          DEALLOCATE(pAmbient%ptr)
00124          NULLIFY(pAmbient%ptr)
00125          CALL ObjectListRemove(id)
00126       ENDIF
00127       NULLIFY(Ambient)
00128    END SUBROUTINE DestroyAmbient
00129 
00130 
00131    SUBROUTINE AmbientGridInit(Info, Ambient)
00132      TYPE(InfoDef) :: Info
00133      TYPE(AmbientDef), POINTER :: Ambient
00134      INTEGER, POINTER, DIMENSION(:,:,:) :: mSs
00135      REAL(KIND=qPREC), POINTER, DIMENSION(:,:) :: offsets
00136      INTEGER :: nOverlaps
00137      CALL CalcPhysicalOverlaps(Info, Ambient%xBounds, mSs, nOverlaps, offsets, IEVERYWHERE, lHydroPeriodic)
00138      IF (nOverlaps > 0) THEN
00139         CALL PlaceAmbient(Info, Ambient, nOverlaps, mSs, offsets)
00140         DEALLOCATE(mSs, offsets)
00141      END IF
00142    END SUBROUTINE AmbientGridInit
00143 
00144 
00145    SUBROUTINE AmbientBeforeStep(Info, Ambient)
00146      TYPE(InfoDef) :: Info
00147      TYPE(AmbientDef), POINTER :: Ambient
00148      INTEGER, POINTER, DIMENSION(:,:,:) :: mSs
00149      REAL(KIND=qPREC), POINTER, DIMENSION(:,:) :: offsets
00150      INTEGER :: nOverlaps
00151      INTEGER :: i,j
00152      DO i=1,nDim
00153         DO j=1,2
00154            IF (Ambient%PersistInBoundaries(i,j)) THEN
00155               CALL CalcPhysicalOverlaps(Info, Ambient%xBounds, mSs, nOverlaps, offsets, IBOUNDARY(i,j), lHydroPeriodic)
00156               IF (nOverlaps > 0) THEN
00157                  CALL PlaceAmbient(Info, Ambient, nOverlaps, mSs, offsets)
00158                  DEALLOCATE(mSs, offsets)
00159               END IF
00160            END IF
00161         END DO
00162      END DO
00163    END SUBROUTINE AmbientBeforeStep
00164 
00165    SUBROUTINE AmbientSetErrFlag(Info, Ambient)
00166       TYPE(InfoDef) :: Info
00167       Type(AmbientDef), POINTER :: Ambient
00168    END SUBROUTINE AmbientSetErrFlag
00169   
00170 
00173    SUBROUTINE PlaceAmbient(Info, Ambient, nOverlaps, mSs, offsets)
00174      TYPE(InfoDef) :: Info
00175      TYPE(AmbientDef), POINTER :: Ambient
00176      INTEGER :: i, j, k, l, n
00177      REAL(KIND=qPREC) :: pos(3), rho, vel(3), x, y, z, press, dx, r, phi
00178      INTEGER :: location, nOverlaps
00179      INTEGER, POINTER, DIMENSION(:,:,:) :: mSs
00180      INTEGER, DIMENSION(3,2) :: mS
00181      REAL(KIND=qPREC), DIMENSION(:,:), POINTER :: offsets
00182      REAL(KIND=qPREC), DIMENSION(3) :: offset
00183      REAL(KIND=qPREC), DIMENSION(3,2) :: universe
00184 
00185      dx=levels(Info%level)%dx
00186      IF (nOverlaps > 0) THEN
00187         DO n=1, nOverlaps
00188            mS=mSs(n,:,:)
00189            offset=offsets(n,:)
00190            IF(MaintainAuxArrays) THEN
00191               DO i=1,nAux
00192                  mS(i,2)=mS(i,2)+1
00193                  Info%aux(mS(1,1):mS(1,2),mS(2,1):mS(2,2),mS(3,1):mS(3,2),i)=Ambient%B(i)
00194                  mS(i,2)=mS(i,2)-1
00195               END DO
00196            END IF
00197 
00198            DO i=mS(1,1), mS(1,2)
00199               x=Info%xBounds(1,1)+(real(i)-half)*dx + offset(1)
00200               DO j=mS(2,1), mS(2,2)
00201                  y=Info%xBounds(2,1)+(real(j)-half)*dx + offset(2)
00202                  DO k=mS(3,1), mS(3,2)
00203                     z=Info%xBounds(3,1)+(real(k)-half)*merge(dx, 0d0, nDim==3) + offset(3)
00204                     pos=(/x,y,z/)
00205                     rho=Ambient%Density
00206                     press=Ambient%pressure
00207                     vel=Ambient%velocity
00208                     IF (Associated(Ambient%DensityPerturbation)) rho = Ambient%Density+PerturbationValue(Ambient%DensityPerturbation,pos)
00209                     IF (ASSOCIATED(Ambient%VelocityPerturbation)) vel = Ambient%velocity+VectorPerturbationValue(Ambient%VelocityPerturbation,pos)
00210                     IF (ASSOCIATED(Ambient%PointGravity)) THEN
00211                        r = sqrt(sum(pos(1:nDim)**2))
00212                        phi=GravityPotential(Ambient%PointGravity%mass, pos, Ambient%PointGravity%soft_length, Ambient%PointGravity%soft_function)
00213                        rho=min(exp(-phi/(ambient%pressure/ambient%density)+log(ambient%density)), 1d300)
00214                        press=rho*ambient%pressure/ambient%density
00215 
00216                        !                        write(*,*) phi, ambient%pressure/ambient%density, phi/(ambient%pressure/ambient%density), exp(phi/(ambient%pressure/ambient%density))
00217                        !                        STOP
00218 
00219                     ELSEIF (ASSOCIATED(Ambient%Profile)) THEN
00220                        r = sqrt(sum(pos(1:nDim)**2))
00221                        DO l=1, size(Ambient%Profile%fields)
00222                           IF (Ambient%Profile%fields(l) == Mass_Field) rho=GetProfileValue(r,Mass_Field,Ambient%Profile)
00223                           IF (Ambient%Profile%fields(l) == P_Field) press=GetProfileValue(r,P_Field,Ambient%Profile)
00224                        END DO
00225 
00226                     END IF
00227                     Info%q(i,j,k,1:NrHydroVars)=0
00228                     Info%q(i,j,k,1)=rho
00229                     IF (ivx /= 0) Info%q(i,j,k,ivx)=vel(1)
00230                     IF (ivy /= 0) Info%q(i,j,k,ivy)=vel(2)
00231                     IF (ivz /= 0) Info%q(i,j,k,ivz)=vel(3)
00232                     IF (iE /= 0) Info%q(i,j,k,iE)=press
00233                     IF (iBx /= 0) Info%q(i,j,k,iBx)=Ambient%B(1)
00234                     IF (iBy /= 0) Info%q(i,j,k,iBy)=Ambient%B(2)
00235                     IF (iBz /= 0) Info%q(i,j,k,iBz)=Ambient%B(3)     
00236                     IF (Ambient%iTracer > 0) Info%q(i,j,k,Ambient%iTracer) = Ambient%density
00237                  END DO
00238               END DO
00239            END DO
00240            CALL Prim_to_cons(Info, mS)
00241         END DO
00242      END IF
00243    END SUBROUTINE PlaceAmbient
00244 
00245  END MODULE Ambients
00246 
 All Classes Files Functions Variables