Scrambler
1
|
00001 !######################################################################### 00002 ! 00003 ! Copyright (C) 2003-2012 Department of Physics and Astronomy, 00004 ! University of Rochester, 00005 ! Rochester, NY 00006 ! 00007 ! physics_control.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 00028 00031 00035 00038 MODULE PhysicsControl 00039 00040 USE GlobalDeclarations 00041 USE EllipticDeclarations 00042 USE PhysicsDeclarations 00043 USE SlopeLim 00044 USE SourceControl 00045 00046 IMPLICIT NONE 00047 PRIVATE 00048 PUBLIC PhysicsInit, PhysicsFinalizeTracers, PhysicsFinalizeInit 00049 00050 CONTAINS 00051 00053 SUBROUTINE PhysicsInit() 00054 00055 INTEGER i,iErr 00056 INTEGER :: Correctmbc(3,2) 00057 ! Open and read in data from the physics data file. 00058 OPEN(UNIT=PHYSICS_DATA_HANDLE, FILE=PHYSICS_DATA_FILE, STATUS='old', FORM='formatted') 00059 READ(PHYSICS_DATA_HANDLE, NML=PhysicsData) 00060 00061 ! Check the correctness of MANDATORY and OPTIONAL parameter initialization 00062 ! [BDS][20100707]: Commented out source-term version. 00063 00064 IF (gamma <= zero) THEN 00065 IF (MPI_ID == 0) write(*,*) 'Gamma must be larger than zero' 00066 STOP 00067 END IF 00068 00069 ! Initialize gamma-related constants for the exact Riemann solver 00070 gamma1 = gamma-1.0 00071 gamma2 = 2.0/gamma1 00072 gamma3 = 0.5*gamma1/gamma 00073 gamma4 = (gamma+1.)/gamma1 00074 gamma5 = 2.0/SQRT(2.0*gamma*gamma1) 00075 gamma6 = 1.0/gamma 00076 gamma7 = 1.0/gamma1 00077 gamma8 = (gamma+1)/(2.0*gamma) 00078 gamma9 = 2.0/(gamma+1.0) 00079 gamma10 = 2.0*gamma/gamma1 00080 gamma11 = gamma1/(gamma+1.0) 00081 gamma12 = 0.5*gamma1 00082 gamma13 = 2.0-gamma 00083 gamma14 = gamma1/gamma 00084 gamma15 = gamma/gamma1 00085 00086 !Currently no support for source terms. 00087 IF(nDim > 2 .AND. iCylindrical /= NoCyl) THEN 00088 IF (MPI_ID.eq.0) PRINT *,'!!! physics_control ERROR: Cylindrical symmetry,', & 00089 'is not implemented for nDim > 2' 00090 STOP 00091 END IF 00092 00093 IF (iCylindrical /= NoCyl) THEN 00094 IF (lMHD) THEN 00095 IF (iCylindrical /= WithAngMom) THEN 00096 IF (MPI_ID.eq.0) PRINT *, 00097 '!!! setprob WARNING: Setting iCylindrical from & NoAngMom to WithMom since lMHD == T' 00098 iCylindrical = WithAngMom 00099 END IF 00100 END IF 00101 IF (ANY(GxBounds(1,:).eq.(/0.0,0.0/) .AND. Gmthbc(1,:) /= (/REFLECT_CYLINDRICAL,REFLECT_CYLINDRICAL/))) THEN 00102 Correctmbc=Gmthbc 00103 WHERE(GxBounds(1,:) == 0.) CorrectMbc(1,:)=REFLECT_CYLINDRICAL 00104 IF (MPI_ID.eq.0) print*,'Left or right boundary corresponds to the z axis. So the boundary conditions Gmthbc is being set to', Correctmbc, 'instead of ', Gmthbc 00105 Gmthbc=Correctmbc 00106 END IF 00107 IF (lMHD) THEN 00108 IF (ANY(GxBounds(2,:).eq.0. .AND. Gmthbc(2,:) /= REFLECT_BPARALLEL)) THEN 00109 Correctmbc=Gmthbc 00110 WHERE(GxBounds(2,:) == 0) CorrectMbc(2,:)=REFLECT_BPARALLEL 00111 IF (MPI_ID.eq.0) print*,'Top or bottom boundary corresponds to the r axis. Perhaps they should be', Correctmbc, 'instead of ', Gmthbc 00112 END IF 00113 ELSE 00114 IF (ANY(GxBounds(2,:).eq.0. .AND. (Gmthbc(2,:) /= REFLECT_BPARALLEL .OR. Gmthbc(2,:) /= REFLECT_WALL))) THEN 00115 Correctmbc=Gmthbc 00116 WHERE(GxBounds(2,:) == 0) CorrectMbc(2,:)=REFLECT_BPARALLEL 00117 IF (MPI_ID.eq.0) print*,'Top or bottom boundary corresponds to the r axis. Perhaps they should be', Correctmbc, 'instead of ', Gmthbc 00118 END IF 00119 END IF!lMHD 00120 END IF!icyl 00121 00122 R2DEff=(sqrt(product(GxBounds(1:2,2)-GxBounds(1:2,1)))/sqrt(pi)*exp(-.5d0)) 00123 ! Source Terms are not implemented in this iteration 00124 00125 CALL SetConservedVariableIndices() 00126 00127 IF (lSelfGravity) THEN 00128 IF (iThreaded /= NON_THREADED) THEN 00129 IF (MPI_ID == 0) THEN 00130 PRINT*, 'SelfGravity only supported with iThreading = -1' 00131 PRINT*, 'Setting iThreading to -1' 00132 END IF 00133 iThreaded = NON_THREADED 00134 LevelBalance=(/0d0,0d0/) 00135 END IF 00136 lSinkParticles=.true. 00137 lParticles=.true. 00138 lElliptic=.true. 00139 END IF 00140 lUniformGravity = UniformGravity /= 0d0 00141 00142 lExplicit = lResistive.or.lConductive.or.lViscous 00143 00144 CALL SetPhysicalScales() 00145 00146 lCanAddTracers=.true. 00147 END SUBROUTINE PhysicsInit 00148 00150 SUBROUTINE PhysicsFinalizeTracers() 00151 !Build Tracers onto q 00152 CALL SrcInitTracers 00153 IF (NrTracerVars > 0) THEN 00154 nTracerHI = NrCons + NrTracerVars 00155 END IF 00156 NrHydroVars=NrCons+NrTracerVars 00157 NrEllipticVars=0 00158 lCanAddTracers=.false. 00159 END SUBROUTINE PhysicsFinalizeTracers 00160 00162 SUBROUTINE PhysicsFinalizeInit() 00163 CALL PhysicsFinalizeTracers() 00164 CALL FinalizeVariableIndices 00165 END SUBROUTINE PhysicsFinalizeInit 00166 00167 00169 SUBROUTINE SetConservedVariableIndices 00170 INTEGER :: i, iEllipticFree, iDiagnostics 00171 00172 CALL AddFieldVar(irho, "rho") 00173 CALL AddFieldVar(ivx, "px") 00174 IF (nDim >= 2 .OR. lMHD) CALL AddFieldVar(ivy, "py") 00175 IF (nDim == 3 .OR. lMHD .OR. iCylindrical == WithAngMom) CALL AddFieldVar(ivz, "pz") 00176 IF (iEOS /= EOS_ISOTHERMAL) CALL AddFieldVar(iE, "E") 00177 IF (lMHD) THEN 00178 CALL AddFieldVar(iBx, "Bx") 00179 CALL AddFieldVar(iBy, "By") 00180 CALL AddFieldVar(iBz, "Bz") 00181 END IF 00182 mydim = merge(nDim, 3, iCylindrical == NOCYL) 00183 GravityDim = mydim 00184 NrCons=NrFieldVars 00185 m_low=ivx 00186 IF (ivz /= 0) THEN 00187 m_high=ivz 00188 ELSEIF (ivy /= 0) THEN 00189 m_high=ivy 00190 ELSE 00191 m_high=ivx 00192 END IF 00193 00194 iAngMom=ivz 00195 imom(1:3)=(/ivx,ivy,ivz/) 00196 iB(1:3)=(/iBx,iBy,iBz/) 00197 00198 nTracerLO = NrCons+1 00199 nTracerHI = 0 00200 NrTracerVars = 0 00201 00202 IF (lMHD .AND. nDim >= 2) THEN 00203 nAux=nDim 00204 IF (nDim == 2) THEN 00205 nEmf=1 00206 iEx=0 00207 iEy=0 00208 iEz=1 00209 EmfLoc=(/0,0,1/) 00210 EmfDir=(/3,0,0/) 00211 ELSE 00212 nEmf=3 00213 iEx=1 00214 iEy=2 00215 iEz=3 00216 EmfLoc=(/1,2,3/) 00217 EmfDir=(/1,2,3/) 00218 END IF 00219 MaintainAuxArrays=.true. 00220 ELSE 00221 nAux=0 00222 END IF 00223 00224 END SUBROUTINE SetConservedVariableIndices 00225 00226 00227 SUBROUTINE FinalizeVariableIndices() 00228 USE EllipticDeclarations 00229 INTEGER :: i 00230 ! include elliptics from source terms 00231 ! CALL SrcInitElliptics(NrVars,NrEllipticVars) 00232 nEllipticLo=NrHydroVars+1 00233 nEllipticHi=0 00234 IF (lSelfGravity) THEN 00235 CALL AddElliptic(iPhiGas, 'Gas Phi') 00236 CALL AddElliptic(iPhiDot, 'PhiDot') 00237 IF (iE /= 0) lStoreMassFlux=.true. 00238 END IF 00239 NrPhysicalVars=NrHydroVars+NrEllipticVars 00240 00241 IF (lCheckDivergence) THEN 00242 IF (.NOT. lMHD) THEN 00243 IF (MPI_ID == 0) write(*,*) 'MHD not true in physics.data so ignoring CheckDivergence' 00244 lCheckDivergence=.false. 00245 ELSE 00246 END IF 00247 END IF 00248 00249 !Build BackupVars 00250 NrBackupVars=0 00251 IF (lSelfGravity) THEN 00252 ! irhoOld=NrPhyiscalVars+NrBackupVars+1 00253 ! ivOld(1:nDim)=(/(NrPhysicalVars+NrBackupVars+1+i,i=1,nDim)/) 00254 ! NrBackupVars=NrBackupVars+nDim+1 00255 END IF 00256 NrVars=NrPhysicalVars+NrBackupVars 00257 00258 ! Stuff for the AMR Engine 00259 nEllipticTransfers=0 00260 elliptic_mbc=0 00261 IF (lSelfGravity) THEN 00262 ! nEllipticTransfers=NrEllipticVars !For now transfer all elliptic variables 00263 ! ALLOCATE(EllipticTransferFields(nEllipticTransfers)) 00264 ! EllipticTransferFields=(/(nEllipticLo+i-1, i=1, NrEllipticVars)/) 00265 ! elliptic_mbc=1 !Need one more ghost zone so we can take gradient of phi 00266 EGVars=1!NrEllipticVars 00267 ALLOCATE(EGCopyFields(EGVars)) 00268 EGCopyFields=(/iPhiGas/)!(/(GVars+i,i=1,EGVars)/) 00269 TDVars=1 00270 ALLOCATE(TimeDerivFields(1)) 00271 TimeDerivFields=(/iPhiDot/) 00272 ELSE 00273 EGVars=0 00274 END IF 00275 00276 afterstep_mbc=0 00277 IF(lResistive)THEN 00278 afterstep_mbc=afterstep_mbc+1 00279 END IF 00280 IF(lConductive)THEN 00281 afterstep_mbc=afterstep_mbc+1 00282 END IF 00283 IF(lViscous)THEN 00284 afterstep_mbc=afterstep_mbc+1 00285 END IF 00286 00287 GVars=NrHydroVars 00288 ALLOCATE(GCopyFields(GVars)) 00289 GCopyFields=(/(i,i=1,GVars)/) 00290 00291 00292 ! EGVars=0 00293 ! EGCopyF 00294 ! CopyFields=NrHydroVars+NrEllipticVars 00295 ! ALLOCATE(HydroCopyFields(NrHydroVars)) 00296 00297 ! CopyFields(1:nCopyFields)=(/(i,i=1,nCopyFields)/) 00298 IF (nAux > 0) THEN 00299 nFlux=NrHydroVars-nAux 00300 nProlongate=NrHydroVars+NrEllipticVars-nAux 00301 nRestrict=NrHydroVars-nAux 00302 ALLOCATE(AuxFields(nAux)) 00303 ALLOCATE(ProlongateFields(nProlongate)) 00304 ALLOCATE(RestrictFields(nRestrict)) 00305 ALLOCATE(FluxFields(nFlux)) 00306 AuxFields=(/(i,i=iBx,iBx-1+nAux)/) 00307 ProlongateFields(1:NrPhysicalVars-nAux)=(/(i,i=1,AuxFields(1)-1),(i,i=AuxFields(nAux)+1,NrPhysicalVars)/) !Removes auxfields 00308 RestrictFields(1:NrHydroVars-nAux)=(/(i,i=1,AuxFields(1)-1),(i,i=AuxFields(nAux)+1,NrHydroVars)/) !Removes auxfields 00309 FluxFields(1:NrHydroVars-nAux)=(/(i,i=1,AuxFields(1)-1),(i,i=AuxFields(nAux)+1,NrHydroVars)/) !Removes auxfields 00310 ELSE 00311 nFlux=NrHydroVars 00312 nProlongate=NrPhysicalvars 00313 nRestrict=NrHydroVars 00314 ALLOCATE(ProlongateFields(nProlongate)) 00315 ALLOCATE(RestrictFields(nRestrict)) 00316 ALLOCATE(FluxFields(nFlux)) 00317 ProlongateFields=(/(i,i=1,NrPhysicalVars)/) 00318 RestrictFields=(/(i,i=1,NrHydroVars)/) 00319 FluxFields=(/(i,i=1,NrHydroVars)/) 00320 END IF 00321 00322 !Setup reverse mapping of flux fields 00323 ALLOCATE(invFluxFields(1:NrHydroVars)) 00324 invFluxFields=0 00325 DO i=1,nFlux 00326 invFluxFields(FluxFields(i))=i 00327 END DO 00328 00329 ALLOCATE(InterpMethod(NrPhysicalVars)) 00330 InterpMethod=InterpOpts(1:NrPhysicalVars) 00331 IF (MPI_ID == 0) THEN 00332 IF (nAux > 0) THEN 00333 IF (ANY(InterpMethod(AuxFields(1:nDim)) /= 0)) THEN 00334 InterpMethod(AuxFields(1:nDim))=0 00335 write(*,*) 'Warning: InterpMethod should be set to 0 for aux fields to guarantee truncation limits on divergence' 00336 write(*,'(A,20I3)') 'Consider setting InterpOpts in physics.data to ', InterpMethod 00337 InterpMethod=InterpOpts(1:NrPhysicalVars) 00338 END IF 00339 END IF 00340 END IF 00341 IF (iPhiGas > 0) InterpMethod(iPhiGas)=PARABOLIC_INTERP 00342 IF (iPhiDot > 0) InterpMethod(iPhiDot)=PARABOLIC_INTERP 00343 00344 END SUBROUTINE FinalizeVariableIndices 00345 00347 SUBROUTINE SetPhysicalScales 00348 00349 IF (nScale <= zero .AND. rScale <= zero) THEN 00350 IF (MPI_ID == 0) PRINT*, 'nScale or rScale must be greater than zero' 00351 STOP 00352 END IF 00353 IF (TempScale <= zero .AND. pScale <= zero) THEN 00354 IF (MPI_ID == 0) PRINT*, 'TempScale or pScale must be greater than zero' 00355 STOP 00356 END IF 00357 IF (lScale <= zero) THEN 00358 IF (MPI_ID == 0) PRINT*, 'lScale must be greater than zero' 00359 STOP 00360 END IF 00361 IF (nScale > zero .AND. rScale > zero) nScale=0d0 !revert to using rScale 00362 IF (TempScale > zero .AND. pScale > zero) TempScale = 0d0 !revert to using pScale 00363 IF (Xmu == 0d0 .AND. nScale == 0d0 .AND. TempScale == 0d0) THEN ! Dimensionless run. 00364 TempScale=1d0 00365 nScale=1d0 !This shouldn't be used if Xmu = 0 00366 ELSE 00367 IF (nScale == zero) nScale = rScale/(hMass * Xmu) 00368 IF (rScale == zero) rScale = nScale * hMass * Xmu 00369 IF (pScale == zero) pScale = ((rScale/hMass)/Xmu)* TempScale * Boltzmann 00370 IF (TempScale == zero) TempScale = pScale*Xmu*hMass/(rScale*Boltzmann) !TempScale 00371 END IF 00372 00373 VelScale = SQRT(pScale/rScale) 00374 TimeScale = lScale/VelScale 00375 BScale=SQRT(4d0 * Pi * pScale) 00376 mScale=rScale*lScale**3 00377 ScaleGrav=G*rscale*TimeScale**2 !computational value of 'G' 00378 00379 MinTemp=MinTemp/TempScale !Converts MinTemp to computational units 00380 00381 IF (MinTemp <= 0) THEN 00382 IF (MPI_ID == 0) PRINT*, 'MinTemp must be positive in physics.data. Setting to 1d-10' 00383 MinTemp=1d-10 00384 END IF 00385 IF (MinDensity <= 0) THEN 00386 IF (MPI_ID == 0) PRINT*, 'MinDensity must be positive in physics.data. Setting to 1d-10' 00387 MinDensity=1d-10 00388 END IF 00389 Iso_Speed2=gamma * MinTemp 00390 Iso_Speed=sqrt(Iso_Speed2) 00391 IF (MPI_ID == 0) THEN 00392 open(UNIT=SCALES_DATA_HANDLE, file=SCALES_DATA_FILE, status="unknown") 00393 WRITE(SCALES_DATA_HANDLE, NML=ScalesData) 00394 CLOSE(SCALES_DATA_HANDLE) 00395 END IF 00396 END SUBROUTINE SetPhysicalScales 00397 00398 END MODULE PhysicsControl 00399