Scrambler
1
|
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 GravoTurbulence 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 00036 MODULE Problem 00037 00038 USE DataDeclarations 00039 USE GlobalDeclarations 00040 USE PhysicsDeclarations 00041 USE Ambients 00042 USE Perturbation 00043 USE CoolingSrc 00044 USE Clumps 00045 USE Totals 00046 USE PDFs 00047 USE Projections 00048 USE Fields 00049 USE Histograms 00050 USE ParticleDeclarations 00051 USE Shapes 00052 USE EOS 00053 USE LayoutDeclarations 00054 USE LayoutControl 00055 USE ObjectControl 00056 USE ProcessingDeclarations 00057 IMPLICIT NONE 00058 SAVE 00059 PRIVATE 00060 00061 PUBLIC ProblemModuleInit, ProblemGridInit, & 00062 ProblemBeforeStep, ProblemAfterStep, ProblemSetErrFlag, ProblemBeforeGlobalStep 00063 REAL(KIND=qPREC), DIMENSION(0:MaxDepth) :: cells_per_cooling_length=0 00064 TYPE(CoolingDef),POINTER :: coolingobj 00065 REAL(KIND=qPrec), PARAMETER :: Y=2d-26 00066 ! REAL(KIND=qPREC), DIMENSION(3), PARAMETER :: offset=(/4d0,-4d0,0d0/)/sqrt(2d0) 00067 TYPE(ClumpDef), POINTER :: clump 00068 TYPE(AmbientDef), POINTER :: Ambient 00069 TYPE(LayoutDef), POINTER :: Layout 00070 REAL(KIND=qPREC), DIMENSION(:,:,:,:), POINTER :: data 00071 INTEGER :: RemapLevel 00072 CONTAINS 00073 00075 SUBROUTINE ProblemModuleInit 00076 REAL(KIND=qPREC) :: rho, p 00077 INTEGER :: i, nDwaves, nVwaves, kmin, kmax, nclumps 00078 REAL(KIND=qPREC) :: wavevector(3), amplitude, phase, amplitudes(3), phases(3), JeansLength, temp, density_amp, velocity_amp, beta,chi,pos(3),radius, pmass, scale, t_ff, t_c, velocity_factor, mass_factor, GE, KE, TE, mass_clump, alpha_virial=0d0, thickness 00079 LOGICAL :: lCooling 00080 TYPE(particledef), POINTER :: particle 00081 TYPE(HistogramDef), POINTER :: HISTOGRAM 00082 TYPE(PDFDef), POINTER :: PDF 00083 TYPE(ProjectionDef), POINTER :: Projection 00084 INTEGER :: nbins=100 00085 INTEGER :: nParticleCells=0 00086 NAMELIST /ProblemData/ rho, temp, lCooling, kmin,kmax,beta, density_amp, velocity_amp, cells_per_cooling_length,chi,pos,radius, pmass, scale, nParticleCells, nbins, velocity_factor, mass_factor, alpha_virial, CellsPerJeansLength, RemapLevel, thickness 00087 NAMELIST /WaveData/ wavevector, amplitude, phase 00088 NAMELIST /vWaveData/ wavevector, amplitudes, phases 00089 00090 00091 ! CALL AddTotalVar(TotalMass) 00092 ! CALL AddTotalVar(1, 'mymass') 00093 00094 CALL AddAllTotals(GASCOMP) 00095 CALL AddAllTotals(PARTICLECOMP) 00096 CALL AddAllTotals(BOTHCOMP) 00097 00098 CALL AddDiagnosticVar(MPI_ID_FIELD) 00099 CALL AddDiagnosticVar(ErrFlag_FIELD) 00100 00101 OPEN(UNIT=PROBLEM_DATA_HANDLE, FILE='problem.data') 00102 READ(PROBLEM_DATA_HANDLE,NML=ProblemData) 00103 00104 IF (nParticleCells > 0) THEN 00105 DefaultParticleRadius=nParticleCells 00106 END IF 00107 00108 CALL CreateHistogram(Histogram) 00109 Histogram%Field%iD=1 00110 ! Histogram%Field%name='' 00111 Histogram%Field%component=GASCOMP 00112 Histogram%minvalue=.1d0 00113 Histogram%maxvalue=1d8 00114 Histogram%nbins=nbins 00115 Histogram%scale=LOGSCALE 00116 00117 CALL CreatePDF(PDF) 00118 PDF%Field(1)%iD=Mass_Field 00119 PDF%Field(1)%name='density' 00120 PDF%Field(1)%component=GASCOMP 00121 PDF%Field(2)%iD=P_Field 00122 PDF%Field(2)%name='pressure' 00123 PDF%Field(2)%component=GASCOMP 00124 PDF%minvalue=(/.01d0,100d0/) 00125 PDF%maxvalue=(/1d7,1d6/) 00126 PDF%nbins=(/400,400/) 00127 PDF%Scale=(/LOGSCALE,LOGSCALE/) 00128 PDF%WeightField=VOLUME 00129 00130 CALL CreatePDF(PDF) 00131 PDF%Field(:)%iD=(/Mass_Field, P_Field/) 00132 PDF%Field(1)%name='density' 00133 PDF%Field(2)%name='pressure' 00134 PDF%Field(:)%component=GASCOMP 00135 ! PDF%Field(2)%name='pressure' 00136 PDF%minvalue=(/.01d0,100d0/) 00137 PDF%maxvalue=(/1d7,1d6/) 00138 PDF%nbins=(/400,400/) 00139 PDF%Scale=(/LOGSCALE,LOGSCALE/) 00140 PDF%WeightField=MASS 00141 00142 00143 CALL CreateProjection(projection) 00144 Projection%Field%iD=Mass_Field 00145 Projection%Field%component=BOTHCOMP 00146 Projection%dim=3 00147 ! ALLOCATE(Projection%Image) 00148 ! Projection%Image%Scaling=LOGSCALE 00149 00150 CALL CreateProjection(projection) 00151 Projection%Field%iD=Mass_Field 00152 Projection%Field%component=BOTHCOMP 00153 Projection%dim=1 00154 ! ALLOCATE(Projection%Image) 00155 ! Projection%Image%Scaling=LOGSCALE 00156 00157 CALL CreateProjection(projection) 00158 Projection%Field%iD=Mass_Field 00159 Projection%Field%component=BOTHCOMP 00160 Projection%pow=1d0 00161 Projection%dim=2 00162 ! ALLOCATE(Projection%Image) 00163 ! Projection%Image%Scaling=LOGSCALE 00164 00165 !! CALL CreateProjection(projection) 00166 ! Projection%Field%iD=CoolingStrength_Field 00167 ! Projection%Field%component=GASCOMP 00168 ! Projection%pow=1d0 00169 ! Projection%dim=3 00170 ! ALLOCATE(Projection%Image) 00171 ! Projection%Image%Scaling=LOGSCALE 00172 00173 ! CALL CreateProjection(projection) 00174 ! Projection%Field%iD=CoolingStrength_Field 00175 ! Projection%Field%component=GASCOMP 00176 ! Projection%pow=1d0 00177 ! Projection%dim=2 00178 ! ALLOCATE(Projection%Image) 00179 ! Projection%Image%Scaling=LOGSCALE 00180 00181 ! CALL CreateProjection(projection) 00182 ! Projection%Field%iD=CoolingStrength_Field 00183 ! Projection%Field%component=GASCOMP 00184 ! Projection%pow=1d0 00185 ! Projection%dim=1 00186 ! ALLOCATE(Projection%Image) 00187 ! Projection%Image%Scaling=LOGSCALE 00188 00189 00190 00191 00192 NULLIFY(Ambient) 00193 CALL CreateAmbient(Ambient) 00194 IF (lCooling) THEN 00195 IF (.NOT. lRestart) THEN 00196 CALL CreateCoolingObject(coolingobj) 00197 ELSE 00198 coolingobj => firstcoolingobj 00199 END IF 00200 coolingobj%iCooling=IICOOL 00201 coolingobj%floortemp=1d0 00202 coolingobj%mintemp=0.001 00203 END IF 00204 IF (iE /= 0 .AND. temp == 0) THEN 00205 IF (lCooling) THEN 00206 CALL InitIICool(coolingobj) 00207 temp=GetIICoolEqTemp(rho) 00208 ! Ambient%Pressure=Ambient%density*temp/TempScale 00209 IF (MPI_ID ==0) write(*,*) 'Equilibrium temp (K)=', temp 00210 ELSE 00211 IF (MPI_ID == 0) write(*,*) 'Temperature must be specified if cooling is not used' 00212 STOP 00213 END IF 00214 END IF 00215 IF (iE == 0) temp=Iso_Speed2*TempScale 00216 Ambient%density=rho/chi 00217 IF (iE /= 0) Ambient%pressure=rho*temp/TempScale !Ambient%density*GetIICoolEqTemp(Ambient%density)/TempScale 00218 IF (MPI_ID == 0) THEN 00219 ! IF (iE /= 0) write(*,*) 'Ambient temp=', GetIICoolEqTemp(Ambient%density) 00220 IF (iE /= 0) write(*,*) 'clump temp=', temp 00221 END IF 00222 00223 IF (lSelfGravity) THEN 00224 JeansLength=sqrt(gamma*temp/TempScale*Pi/ScaleGrav/rho) 00225 mass_clump=4d0/3d0*mass_factor*pi*radius**3*rho 00226 GE=3d0/5d0*mass_clump**2*ScaleGrav/radius 00227 TE=mass_clump*temp/TempScale*3/2 00228 IF (alpha_virial /= 0) THEN 00229 KE=alpha_virial*GE/2-TE 00230 IF (KE < 0) THEN 00231 write(*,*) 'virial parameter too low given thermal energy' 00232 STOP 00233 END IF 00234 velocity_amp=sqrt(2d0*KE/mass_clump)/velocity_factor 00235 ! KE=.5*mass_clump*(velocity_amp*velocity_factor)**2 00236 ! alpha_virial=2*(KE+TE)/GE 00237 IF (MPI_ID == 0) write(*,*) 'velocity_amp=', velocity_amp 00238 ELSE 00239 KE=.5*mass_clump*(velocity_amp*velocity_factor)**2 00240 alpha_virial=2*(KE+TE)/GE 00241 00242 END IF 00243 IF (MPI_ID == 0) THEN 00244 write(*,'(A30,2A25)') 'Quantity: ', 'Computationl Units', 'Astro Units' 00245 write(*,'(A30,2E25.16)') 'Temperature of sphere', temp/TempScale, temp 00246 write(*,'(A30,2E25.16)') 'sound speed', sqrt(gamma*temp/TempScale), sqrt(gamma*temp/TempScale)*VelScale/1e5 00247 write(*,'(A30,2E25.16)') 'free fall time', sqrt(3*pi/32/ScaleGrav/rho), sqrt(3*pi/32/ScaleGrav/rho)*TimeScale/yr/1e6 00248 write(*,'(A30,2E25.16)') 'Radius of sphere', radius, radius*lScale/pc 00249 write(*,'(A30,2E25.16)') 'Mass of sphere', mass_clump, mass_clump*mScale/msolar 00250 WRITE(*,'(A30,2E25.16)') 'Jeans Length ', JeansLength, JeansLength*lScale/pc 00251 WRITE(*,'(A30,2E25.16)') 'Jeans Mass ', 4d0/3d0*pi*(JeansLength/2d0)**3*rho, 4d0/3d0*pi*(JeansLength/2d0)**3*rho*mScale/mSolar 00252 WRITE(*,'(A30,2E25.16)') 'Gravitational Energy', GE, GE*pScale*lScale**3 00253 write(*,'(A30,2E25.16)') 'Kinetic Energy of sphere = ', KE, KE*pScale*lScale**3 00254 write(*,'(A30,2E25.16)') 'Thermal Energy of sphere = ', TE, TE*pScale*lScale**3 00255 write(*,'(A30,2E25.16)') 'Velocity dispersion= ', (KE+TE)/mass_clump, (KE+TE)/mass_clump*velScale/1e5 00256 write(*,*) 'Virial parameter = ', 2*(KE+TE)/GE 00257 END IF 00258 00259 END IF 00260 00261 00262 DO i=0,-1 00263 rho=10d0**(.25*REAL(i)-1d0) 00264 IF (iE /= 0) THEN 00265 temp=GetIICoolEqTemp(rho) 00266 ELSE 00267 temp=Iso_Speed2*TempScale 00268 END IF 00269 JeansLength=sqrt(gamma*temp/TempScale*Pi/ScaleGrav/rho) !pc 00270 t_ff=sqrt(3d0*pi/32d0/ScaleGrav/rho)*TimeScale/yr/1e6 !myr 00271 t_c=rho*temp/max(abs(rho*nScale*Y * (rho*nScale*IICoolingRate(temp*TempScale) ) * coolingobj%ScaleCool * gamma1 ),1d-6)*TimeScale/yr/1e6 00272 00273 write(*,'(8E24.15)') rho, temp, JeansLength, t_ff, t_c 00274 END DO 00275 ! STOP 00276 IF (pMass > 0) THEN 00277 CALL CreateParticle(particle) 00278 Particle%Q(1)=pMass 00279 Particle%xloc=pos !+offset 00280 Particle%radius=1 00281 CALL AddSinkParticle(particle) 00282 CALL CreatePointGravityObject(particle%PointGravityObj) 00283 particle%PointGravityObj%mass=pMass 00284 particle%PointGravityObj%x0(1:nDim)=particle%xloc 00285 particle%PointGravityObj%soft_length=4d0*sink_dx !particle%radius*sink_dx 00286 particle%PointGravityObj%soft_function=SPLINESOFT 00287 particle%iAccrete=FEDERRATH_ACCRETION 00288 ! NULLIFY(particle) 00289 ! CALL CreateParticle(particle) 00290 ! Particle%Q(1)=pMass 00291 ! Particle%xloc=pos-offset 00292 ! Particle%radius=1 00293 ! CALL AddSinkParticle(particle) 00294 00295 ELSE 00296 CALL CreateClump(clump) 00297 clump%radius=radius 00298 clump%density=rho 00299 ! clump%pressure=rho*temp/TempScale 00300 clump%temperature=temp/TempScale 00301 clump%position=pos!+offset 00302 clump%subsample=0 00303 clump%thickness=thickness 00304 CALL UpdateClump(clump) 00305 ! NULLIFY(Clump) 00306 ! CALL CreateClump(clump) 00307 ! clump%radius=radius 00308 ! clump%density=rho 00309 ! clump%pressure=rho*temp/TempScale 00310 ! clump%temperature=temp/TempScale 00311 ! clump%position=pos-offset 00312 ! CALL UpdateClump(clump) 00313 END IF 00314 00315 ! IF (nDwaves > 0) THEN 00316 ! CALL CreatePerturbationObj(Ambient%DensityPerturbation) 00317 ! Ambient%DensityPerturbation%type=COSINESERIES 00318 ! CALL InitPerturbationWaves(Ambient%DensityPerturbation, nDwaves) 00319 ! DO i=1,nDwaves 00320 ! READ(PROBLEM_DATA_HANDLE, NML=WaveData) 00321 ! CALL AddPerturbationWave(Ambient%DensityPerturbation, wavevector*2d0*Pi/max(GxBounds(:,2)-GxBounds(:,1),1!e-6), phase, amplitude) 00322 ! END DO 00323 ! END IF 00324 ! IF (density_amp > 0d0) THEN 00325 ! CALL CreatePerturbationObj(Clump%DensityPerturbation) 00326 ! CALL CreateSpectra(Clump%DensityPerturbation, kmin, kmax, beta, density_amp, scale) 00327 ! END IF 00328 ! IF (velocity_amp > 0) THEN 00329 ! CALL CreateVectorPerturbationObj(Clump%VelocityPerturbation) 00330 ! CALL CreateSolenoidalSpectra(Clump%VelocityPerturbation, kmin, kmax, beta, velocity_amp, scale) 00331 ! END IF 00332 00333 ! IF (nVwaves > 0) THEN 00334 ! CALL CreateVectorPerturbationObj(Ambient%VelocityPerturbation) 00335 ! Ambient%VelocityPerturbation%type=COSINESERIES 00336 ! CALL InitVectorPerturbationWaves(Ambient%VelocityPerturbation, nVwaves) 00337 ! DO i=1,nVwaves 00338 ! READ(PROBLEM_DATA_HANDLE, NML=VWaveData) 00339 ! CALL AddVectorPerturbationWave(Ambient%VelocityPerturbation, wavevector*2d0*Pi/max(GxBounds(:,2)-GxBounds!(:,1),1e-6), phases, amplitudes) 00340 ! END DO 00341 ! END IF 00342 ! write(*,*) VectorPerturbationValue(Ambient%VelocityPerturbation,(/.1d0,.1d0,.1d0/)) 00343 00344 CLOSE(PROBLEM_DATA_HANDLE) 00345 00346 END SUBROUTINE ProblemModuleInit 00347 00350 SUBROUTINE ProblemGridInit(Info) 00351 !! @brief Initializes the grid data according to the requirements of the problem. 00352 !! @param Info A grid structure. 00353 TYPE (InfoDef) :: Info 00354 REAL(KIND=qPREC) :: p, f, rho, vel(3), pos(3) 00355 INTEGER :: i,j,k 00356 IF (levels(Info%level)%tnow == start_Time .AND. lRestart) THEN 00357 IF (Info%level < RemapLevel) THEN 00358 ! CALL AmbientGridInit(Info) 00359 ! CALL ClumpGridInit(Info) 00360 CALL ObjectsGridInit(Info,AMBIENTOBJ) 00361 CALL ObjectsGridInit(Info,CLUMPOBJ) 00362 ELSE 00363 DO i=1,Info%mX(1) 00364 DO j=1,Info%mx(2) 00365 DO k=1,Info%mx(3) 00366 IF (IsInShape(Clump%Shape, CellPos(Info, i, j, k), pos)) THEN 00367 f=smooth_tanh(sqrt(sum(pos**2))/Clump%Shape%size_param(1),Clump%thickness) 00368 p=Press(Info%q(i,j,k,:))*f+Ambient%pressure*(1d0-f) 00369 rho=Info%q(i,j,k,1)*f+Ambient%density*(1d0-f) 00370 vel(1:nDim)=Info%q(i,j,k,imom(1:nDim))/Info%q(i,j,k,1)*f 00371 ELSE 00372 rho=Ambient%density 00373 p=Ambient%pressure 00374 vel(1:nDim)=0d0 00375 END IF 00376 00377 Info%q(i,j,k,1)=rho 00378 Info%q(i,j,k,imom(1:nDim))=rho*vel(1:nDim) 00379 IF (iE /= 0d0) Info%q(i,j,k,iE)=gamma7*p+half*sum(Info%q(i,j,k,imom(1:nDim)**2))/Info%q(i,j,k,1) 00380 END DO 00381 END DO 00382 END DO 00383 END IF 00384 END IF 00385 END SUBROUTINE ProblemGridInit 00386 00389 SUBROUTINE ProblemBeforeStep(Info) 00390 !! @brief Performs any tasks required before the advance step. 00391 !! @param Info A grid structure. 00392 TYPE (InfoDef) :: Info 00393 END SUBROUTINE ProblemBeforeStep 00394 00397 SUBROUTINE ProblemAfterStep(Info) 00398 !! @brief Performs any post-step corrections that are required. 00399 !! @param Info A grid structure. 00400 TYPE (InfoDef) :: Info 00401 END SUBROUTINE ProblemAfterStep 00402 00405 SUBROUTINE ProblemSetErrFlag(Info) 00406 !! @brief Sets error flags according to problem-specific conditions.. 00407 !! @param Info A grid structure. 00408 TYPE (InfoDef) :: Info 00409 REAL(KIND=qPREC) :: P, T, cooling_length, temp, pos(3), dx, dz, duh(3), cool_rate, x, dist_from_boundary 00410 INTEGER :: i,j,k, mO(3,2) 00411 00412 IF (cells_per_cooling_length(Info%level) > 0d0) THEN 00413 temp=levels(Info%level)%dx*cells_per_cooling_length(Info%level) 00414 DO i=1,Info%mX(1) 00415 DO j=1,Info%mX(2) 00416 DO k=1,Info%mX(3) 00417 P=Press(Info%q(i,j,k,:)) 00418 T=TempScale*P/Info%q(i,j,k,1) 00419 cool_rate=max(abs(Info%q(i,j,k,1)*nScale*Y * (Info%q(i,j,k,1)*nScale*IICoolingRate(T) ) * coolingobj%ScaleCool * gamma1 ),1d-6) 00420 cooling_length=sqrt(gamma*P/Info%q(i,j,k,1))*(P/cool_rate) 00421 IF (cooling_length < temp) Info%ErrFlag(i,j,k)=1 00422 END DO 00423 END DO 00424 END DO 00425 END IF 00426 00427 IF (levels(Info%level)%tnow == start_time .AND. lRestart) THEN !remap cloud... 00428 mO=MapBoxToInfo(Clump%Shape%xBounds, Info, 0) 00429 Info%ErrFlag(mO(1,1):mO(1,2), mO(2,1):mO(2,2), mO(3,1):mO(3,2))=1 00430 END IF 00431 END SUBROUTINE ProblemSetErrFlag 00432 00433 SUBROUTINE ProblemBeforeGlobalStep(n) 00434 INTEGER :: n,i 00435 INTEGER :: mB(3,2) 00436 INTEGER, DIMENSION(:), ALLOCATABLE :: FieldID 00437 TYPE(NodeDefList), POINTER :: nodelist 00438 LOGICAL, SAVE :: lFirstTime=.true. 00439 IF (lFirstTime) THEN 00440 IF (iE /= 0) THEN 00441 ALLOCATE(FieldID(5)) 00442 FieldID=(/(i,i=1,5)/) 00443 ELSE 00444 ALLOCATE(FieldID(4)) 00445 FieldID=(/(i,i=1,4)/) 00446 END IF 00447 00448 IF (levels(n)%tnow == start_Time .AND. lRestart .AND. n == 0) THEN 00449 CALL CreateLayout(GmGlobal, layout) 00450 layout%level = 0 00451 mB=layout%mB(MPI_ID,:,:) 00452 ALLOCATE(data(mB(1,1):mB(1,2), mB(2,1):mB(2,2),mB(3,1):mB(3,2),4)) 00453 CALL LoadFieldIntoLayout(layout, data, FieldID) 00454 CALL StoreLayout(layout, data, 'ISO512') 00455 lFirstTime=.false. 00456 ELSEIF (levels(n)%tnow == start_time .AND. lRestart .AND. n == RemapLevel) THEN 00457 layout%level = RemapLevel 00458 layout%mB(:,1:nDim,:)=layout%mB(:,1:nDim,:)+spread(spread(levels(RemapLevel)%mX(1:nDim)/2-GmX(1:nDim)/2,1,MPI_NP), 3, 2) !shift layout to center of level 1 grid 00459 CALL UnloadFieldFromLayout(layout, data, FieldID, lHydroPeriodic, levels(n)%gmbc(1)) 00460 CALL DestroyLayout(layout) 00461 DEALLOCATE(data) 00462 lFirstTime=.false. 00463 END IF 00464 00465 nodelist=>Nodes(n)%p 00466 DO WHILE (ASSOCIATED(nodelist)) 00467 CALL ProblemGridInit(nodelist%self%info) 00468 nodelist=>nodelist%next 00469 END DO 00470 DEALLOCATE(FieldID) 00471 END IF 00472 END SUBROUTINE ProblemBeforeGlobalStep 00473 00474 END MODULE Problem