Scrambler
1
|
00001 !######################################################################### 00002 ! 00003 ! Copyright (C) 2003-2012 Department of Physics and Astronomy, 00004 ! University of Rochester, 00005 ! Rochester, NY 00006 ! 00007 ! module_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 ModuleControl 00039 00040 USE GlobalDeclarations 00041 USE DataDeclarations 00042 USE PhysicsDeclarations 00043 USE Clumps 00044 USE Disks 00045 USE UniformRegions 00046 USE SplitRegions 00047 USE CollidingFlows 00048 USE Winds 00049 USE Ambients 00050 USE Outflows 00051 USE Problem 00052 USE ParticleDeclarations 00053 USE Refinements 00054 USE ObjectControl 00055 IMPLICIT NONE 00056 00057 PUBLIC ModuleObjectsInit, ModuleProblemInit, GridInit, BeforeStep, SetErrFlag, ApplyPhysicalBC, AfterFixup, BeforeGlobalStep 00058 00059 ! NAMELIST /ModulesData/ rhoOut,pOut,vxOut,vyOut,vzOut,BxOut,ByOut,BzOut 00060 CONTAINS 00061 00062 SUBROUTINE BeforeGlobalStep(n) 00063 USE PointGravitySrc 00064 INTEGER :: n 00065 IF (n == MaxLevel) THEN 00066 CALL CheckParticlePositions() 00067 END IF 00068 CALL ClumpBeforeGlobalStep(n) 00069 CALL ProblemBeforeGlobalStep(n) 00070 END SUBROUTINE BeforeGlobalStep 00071 00073 SUBROUTINE ModuleObjectsInit 00074 !CALL InitAmbients 00075 !CALL UpdateClumps 00076 !CALL UpdateDisks 00077 !CALL InitWinds 00078 !CALL InitOutflows 00079 !CALL InitCollidingFlows 00080 !CALL InitUniformRegions 00081 !CALL InitSplitRegions 00082 !CALL InitRefinements 00083 00084 ! Read in modules' data files 00085 END SUBROUTINE ModuleObjectsInit 00086 00087 00089 SUBROUTINE ModuleProblemInit() 00090 INTEGER :: i 00091 !Need to create default refinements based on refinevariablefactor 00092 IF (RefineVariableFactor(1) /= 0d0) CALL AddRefinementCriterion(Mass_Field, RefineVariableFactor(1)) 00093 IF (iE /= 0) THEN 00094 IF (RefineVariableFactor(iE) /=0d0) CALL AddRefinementCriterion(P_Field, RefineVariableFactor(iE)) 00095 END IF 00096 DO i=m_low, m_high 00097 IF (RefineVariableFactor(i) /= 0d0) CALL AddRefinementCriterion(i, RefineVariableFactor(i), RhoSoundSpeed_Field) 00098 END DO 00099 IF (lMHD) THEN 00100 DO i=iBx, iBz 00101 IF (RefineVariableFactor(i) /= 0d0) CALL AddRefinementCriterion(i, RefineVariableFactor(i), SqrtPress_Field) 00102 END DO 00103 END IF 00104 IF (lSelfGravity) THEN 00105 CALL AddRefinementCriterion(GasPotential_Field, RefineVariableFactor(i), SoundSpeed2_Field, LINEARSCALE) 00106 CALL AddRefinementThreshold(JeansLength_Field, LESSTHAN, (/(4d0*levels(i)%dx,i=0,MaxLevel)/)) 00107 END IF 00108 CALL ProblemModuleInit() 00109 ! CALL InitAmbient 00110 END SUBROUTINE ModuleProblemInit 00111 00112 00115 SUBROUTINE GridInit(Info) 00116 TYPE (InfoDef) :: Info ! Data associated with this grid 00117 INTEGER :: mB(3,2), rmbc 00118 00119 mb(ndim+1:3,:)=1 00120 mb(1:nDim,1)=1 00121 mb(1:nDim,2)=Info%mX(1:nDim) 00122 IF (Info%level == 0) THEN 00123 IF (iPhiDot /= 0) Info%q(:,:,:,iPhiDot)=0 00124 IF (iPhiGas /= 0) Info%q(:,:,:,iPhiGas) = 0 00125 END IF 00126 !CALL AmbientGridInit(Info) 00127 !CALL CollidingFlowGridInit(Info) 00128 !CALL ClumpGridInit(Info) 00129 !CALL DiskGridInit(Info) 00130 !CALL WindGridInit(Info) 00131 !CALL OutflowGridInit(Info) 00132 !CALL UniformRegionGridInit(Info) 00133 !CALL SplitRegionGridInit(Info) 00134 00135 ! CALL ObjectsGridInit(Info, AMBIENTOBJ) 00136 ! CALL ObjectsGridInit(Info, COLLIDINGFLOWOBJ) 00137 ! CALL ObjectsGridInit(Info, CLUMPOBJ) 00138 ! CALL ObjectsGridInit(Info, DISKOBJ) 00139 ! CALL ObjectsGridInit(Info, WINDOBJ) 00140 ! CALL ObjectsGridInit(Info, OUTFLOWOBJ) 00141 ! CALL ObjectsGridInit(Info, UNIFORMREGIONOBJ) 00142 ! CALL ObjectsGridInit(Info, SPLITREGIONOBJ) 00143 00144 CALL ObjectsGridInit(Info) 00145 00146 CALL ProblemGridInit(Info) 00147 CALL Protectq(Info, mb, 'grid init', .TRUE.) 00148 CALL CheckDivergence(Info, .TRUE.) 00149 END SUBROUTINE GridInit 00150 00154 SUBROUTINE BeforeStep(Info) 00155 TYPE(InfoDef) :: Info 00156 INTEGER :: i, mB(3,2), rmbc 00157 IF (irhoOld /= 0) Info%q(:,:,:,irhoOld)=Info%q(:,:,:,irho) 00158 DO i=1,nDim 00159 IF (ivOld(i) /= 0) Info%q(:,:,:,ivOld(i))=Info%q(:,:,:,imom(i))/Info%q(:,:,:,irho) 00160 END DO 00161 mb(ndim+1:3,:)=1 00162 rmbc=levels(Info%level)%gmbc(levels(Info%level)%step) 00163 mb(1:nDim,1)=1-rmbc 00164 mb(1:nDim,2)=Info%mX(1:nDim)+rmbc 00165 00166 00167 ! If the floor temperature is set, then make sure the energy density 00168 ! throughout the grid does not drop below the minimum established by 00169 ! MinTemp * rho + 1/2 (rho*|v|^2 + |B|^2). 00170 00171 ! IF(MinTemp > zero .AND. iE /= 0) THEN 00172 ! IF(lMHD) THEN 00173 ! Info%q(:,:,:,iE) = MAX(Info%q(:,:,:,iE), & 00174 ! MinTemp/EOSConstants*Info%q(:,:,:,1)/(gamma1)+ & 00175 ! half*SUM(Info%q(:,:,:,m_low:m_high)**2,DIM=4)/Info%q(:,:,:,1)+& 00176 ! half*SUM(Info%q(:,:,:,iBx:iBz)**2,DIM=4)) 00177 ! ELSE 00178 ! Info%q(:,:,:,iE) = MAX( Info%q(:,:,:,iE), & 00179 ! MinTemp/EOSConstants*Info%q(:,:,:,1)/(gamma1)+ & 00180 ! half*SUM(Info%q(:,:,:,m_low:m_high)**2,DIM=4)/Info%q(:,:,:,1)) 00181 ! END IF 00182 00183 ! END IF 00184 00185 ! set EMF mask on level 0. 00186 ! IF(lMHD .AND. MaintainAuxArrays .AND. Info%level==0) THEN 00187 ! Info%childemf = UNDEFINED 00188 ! END IF 00189 00190 ! *** NO GENERIC BEFORESTEP *** 00191 !CALL Protectq(Info, mb, 'before step') 00192 !IF (lWinds) CALL WindsBeforeStep(Info) 00193 !IF (lOutflows) CALL OutflowBeforeStep(Info) 00194 !CALL UniformRegionBeforeStep(Info) 00195 !CALL SplitRegionBeforeStep(Info) 00196 !CALL CollidingFlowBeforeStep(Info) 00197 !CALL ClumpBeforeStep(Info) 00198 ! *** NO GENERIC BEFORESTEP *** 00199 00200 CALL Protectq(Info, mb, 'before step') 00201 00202 ! CALL ObjectsBeforeStep(Info, AMBIENTOBJ) 00203 ! CALL ObjectsBeforeStep(Info, WINDOBJ) 00204 ! CALL ObjectsBeforeStep(Info, DISKOBJ) 00205 ! CALL ObjectsBeforeStep(Info, OUTFLOWOBJ) 00206 ! CALL ObjectsBeforeStep(Info, UNIFORMREGIONOBJ) 00207 ! CALL ObjectsBeforeStep(Info, SPLITREGIONOBJ) 00208 ! CALL ObjectsBeforeStep(Info, COLLIDINGFLOWOBJ) 00209 ! CALL ObjectsBeforeStep(Info, CLUMPOBJ) 00210 00211 CALL ObjectsBeforeStep(Info) 00212 00213 CALL ProblemBeforeStep(Info) 00214 IF (lCheckDivergence) CALL CheckDivergence(Info) 00215 00216 END SUBROUTINE BeforeStep 00217 00218 00221 SUBROUTINE AfterStep(Info) 00222 00223 TYPE(InfoDef) :: Info 00224 INTEGER :: i, mB(3,2), rmbc 00225 CALL ProblemAfterStep(Info) 00226 mb(ndim+1:3,:)=1 00227 rmbc=levels(Info%level)%ambc(levels(Info%level)%step) 00228 mb(1:nDim,1)=1-rmbc 00229 mb(1:nDim,2)=Info%mX(1:nDim)+rmbc 00230 CALl Protectq(Info, mb, 'after step') 00231 END SUBROUTINE AfterStep 00232 00235 SUBROUTINE AfterFixup(Info) 00236 TYPE(InfoDef) :: Info 00237 END SUBROUTINE AfterFixup 00238 00239 00242 SUBROUTINE SetErrFlag(Info) 00243 USE ParticleInfoOps 00244 USE ParticleDeclarations 00245 TYPE(InfoDef) :: Info 00246 INTEGER :: rmbc, iErr, n, level 00247 INTEGER, DIMENSION(3) :: mxL, mxH 00248 INTEGER :: i,j,k 00249 REAL(KIND=qPrec) :: x, y, z 00250 REAL(KIND=qPREC), DIMENSION(:,:,:), ALLOCATABLE :: qbuf, cs 00251 level=Info%level 00252 rmbc = levels(level)%gmbc(levels(level)%step) 00253 Info%ErrFlag=0 00254 CALL SinkSetErrFlag(Info) 00255 CALL ObjectsSetErrFlag(Info) 00256 CALL ProblemSetErrFlag(Info) 00257 END SUBROUTINE SetErrFlag 00258 00261 SUBROUTINE ApplyPhysicalBC(Info) 00262 ! Interface declarations 00263 TYPE (InfoDef) :: Info 00264 ! Internal declarations 00265 INTEGER rmbc, level,dim,edge,start,ncells 00266 INTEGER, DIMENSION(3,2):: lGmGlobal, mB 00267 level=Info%level 00268 rmbc=levels(level)%gmbc(levels(level)%step) 00269 lGmGlobal(:,1)=GmGlobal(:,1) 00270 lGmGlobal(:,2)=GmGlobal(:,2)*PRODUCT(levels(0:level-1)%CoarsenRatio) 00271 DO dim=1,nDim 00272 IF (lHydroPeriodic(dim)) CYCLE 00273 DO edge=1,2 00274 IF (edge == 1) THEN 00275 start=lGmGlobal(dim,1)-Info%mGlobal(dim,1) !first cell on left boundary 00276 nCells=start-(1-rmbc)+1 00277 ELSE 00278 start=(lGmGlobal(dim,2)+1)-(Info%mGlobal(dim,1)-1) !first ghost cell on right boundary 00279 nCells=Info%mx(dim)+rmbc-start+1 00280 END IF 00281 IF (nCells > 0) CALL SetGhost(dim,edge,start,nCells,Info,Gmthbc(dim,edge),rmbc) 00282 END DO 00283 END DO 00284 00285 CONTAINS 00286 00287 SUBROUTINE SetGhost(dim, edge, start, nCells, Info, mthbc,rmbc) 00288 INTEGER :: dim, start, nCells, mthbc,edge,dir, level, rmbc,i,j,m 00289 TYPE(InfoDef) :: Info 00290 INTEGER, DIMENSION(3,2) :: ip,iq,ir,is 00291 LOGICAL :: lFirst 00292 INTEGER, DIMENSION(2) :: AuxParFields 00293 REAL(KIND=qPREC) :: aux_psign, aux_nsign 00294 INTEGER, DIMENSION(4) :: ReflectVars 00295 INTEGER :: nReflect 00296 IF (nCells == 0) RETURN 00297 ! write(*,*) dim, edge, start, nCells, mthbc 00298 level=Info%level 00299 dir=(-1)**edge !(edge*2-3) !direction of edge (1=>-1, 2=>1) 00300 ip(nDim+1:3,:)=1 00301 iq(nDim+1:3,:)=1 00302 AuxParFields(1:nDim-1)=modulo(dim+(/(i,i=0,nDim-2)/),nDim)+1 00303 00304 ! Stretch bounds by nCells 00305 DO i=1,nDim 00306 IF (i==dim) CYCLE 00307 ip(i,1)=1-rmbc 00308 ip(i,2)=Info%mX(i)+rmbc 00309 iq(i,:)=ip(i,:) 00310 END DO 00311 SELECT CASE (mthbc) 00312 CASE(EXTRAPOLATED_BOUND) 00313 iq(dim,edge)=start+dir*(nCells-1) 00314 iq(dim,3-edge)=start 00315 ip(dim,:)=start-dir 00316 00317 Info%q(iq(1,1):iq(1,2),iq(2,1):iq(2,2),iq(3,1):iq(3,2),1:NrHydroVars)=& 00318 SPREAD(SUM(Info%q(ip(1,1):ip(1,2),ip(2,1):ip(2,2),ip(3,1):ip(3,2),1:NrHydroVars),dim),dim,nCells) 00319 ! write(*,*) "spreading ", ip, "along ",dim,"into", iq 00320 IF (MaintainAuxArrays) THEN 00321 lFirst=.true. 00322 DO j=1,nDim !Copy parallel auxfields and update b-norm 00323 IF (j == dim) CYCLE 00324 iq(dim,edge)=start+dir*(nCells-1) 00325 iq(dim,3-edge)=start 00326 ip(dim,:)=start-dir 00327 ir=ip 00328 is=iq 00329 ir(dim,:)=ip(dim,:) 00330 is(dim,:)=iq(dim,:) 00331 ir(j,2)=ir(j,2)+1 !stretch aux bounds 00332 is(j,2)=is(j,2)+1 00333 Info%aux(is(1,1):is(1,2),is(2,1):is(2,2),is(3,1):is(3,2),j)=& 00334 SPREAD(SUM(Info%aux(ir(1,1):ir(1,2),ir(2,1):ir(2,2),ir(3,1):ir(3,2),j),dim),dim,nCells) 00335 ! write(*,*) "spreading aux field", j, "at ", ir, "along dim", dim, "into", is 00336 ir(j,2)=ir(j,2)-1 00337 is(j,1)=is(j,1)+1 00338 DO i=1,nCells 00339 iq(dim,:)=start+dir*(i-1)+(edge-1)!start+2-i 00340 ip(dim,:)=start+dir*(i-2)+(edge-1)!start+1-i 00341 ir(dim,:)=start+dir*(i-1) 00342 is(dim,:)=start+dir*(i-1) 00343 IF (lFirst) THEN 00344 Info%aux(iq(1,1):iq(1,2),iq(2,1):iq(2,2),iq(3,1):iq(3,2),dim)=& 00345 Info%aux(ip(1,1):ip(1,2),ip(2,1):ip(2,2),ip(3,1):ip(3,2),dim) 00346 ! write(*,*) "copying aux field", dim, "from" , ip, "into", iq 00347 END IF 00348 ! write(*,*) "and adding...", dim, j, iq, ir, is 00349 Info%aux(iq(1,1):iq(1,2),iq(2,1):iq(2,2),iq(3,1):iq(3,2),dim)= & 00350 Info%aux(iq(1,1):iq(1,2),iq(2,1):iq(2,2),iq(3,1):iq(3,2),dim)+ & 00351 dir*(Info%aux(ir(1,1):ir(1,2),ir(2,1):ir(2,2),ir(3,1):ir(3,2),j)- & 00352 Info%aux(is(1,1):is(1,2),is(2,1):is(2,2),is(3,1):is(3,2),j)) 00353 END DO 00354 is(j,:)=is(j,:)-1 00355 lFirst=.false. 00356 END DO 00357 iq(dim,edge)=start+dir*(nCells-1) 00358 iq(dim,3-edge)=start 00359 ip(dim,:)=iq(dim,:)+1 !tart-dir!iq(dim,:)+1 00360 ! write(*,*) "updating fields", auxfields(dim), "at ", iq, "with average of ", iq, ip 00361 Info%q(iq(1,1):iq(1,2),iq(2,1):iq(2,2),iq(3,1):iq(3,2),AuxFields(dim))=half*(& 00362 Info%aux(iq(1,1):iq(1,2),iq(2,1):iq(2,2),iq(3,1):iq(3,2),dim)+& 00363 Info%aux(ip(1,1):ip(1,2),ip(2,1):ip(2,2),ip(3,1):ip(3,2),dim)) 00364 END IF 00365 CASE(REFLECT_WALL, REFLECT_BPARALLEL, REFLECT_CYLINDRICAL) 00366 aux_psign=1d0 00367 aux_nsign=1d0 00368 SELECT CASE (mthbc) 00369 CASE (REFLECT_WALL) 00370 IF (lMHD) THEN 00371 ReflectVars(1:2)=(/imom(dim), iB(dim)/); nReflect=2; aux_nsign=-1d0 00372 ELSE 00373 ReflectVars(1)=imom(dim); nReflect=1 00374 END IF 00375 CASE(REFLECT_BPARALLEL) 00376 IF (lMHD) THEN 00377 ReflectVars(1:3)=(/imom(dim), iB(modulo(dim+(/0,1/),3)+1)/); nReflect=3; aux_psign=-1d0 00378 ELSE 00379 ReflectVars(1)=imom(dim); nReflect=1 00380 END IF 00381 CASE(REFLECT_CYLINDRICAL) 00382 IF (lMHD) THEN 00383 ReflectVars(1:4)=(/imom(dim), imom(3), iB(dim), iB(3)/); nReflect=4; aux_nsign=-1d0 00384 ELSE 00385 IF (iCylindrical.eq.WithAngMom) THEN 00386 ReflectVars(1:2)=(/imom(dim), imom(3)/); nReflect=2 00387 ELSE 00388 ReflectVars(1:1)=(/imom(dim)/); nReflect=1 00389 END IF 00390 END IF 00391 END SELECT 00392 DO j=1,nCells 00393 iq(dim,:)=start+dir*(j-1) 00394 ip(dim,:)=start-dir*(j) 00395 00396 Info%q(iq(1,1):iq(1,2),iq(2,1):iq(2,2),iq(3,1):iq(3,2),1:NrHydroVars)=& 00397 Info%q(ip(1,1):ip(1,2),ip(2,1):ip(2,2),ip(3,1):ip(3,2),1:NrHydroVars) 00398 Info%q(iq(1,1):iq(1,2),iq(2,1):iq(2,2),iq(3,1):iq(3,2),ReflectVars(1:nReflect)) = & 00399 -Info%q(iq(1,1):iq(1,2),iq(2,1):iq(2,2),iq(3,1):iq(3,2),ReflectVars(1:nReflect)) 00400 00401 IF (MaintainAuxArrays) THEN 00402 DO i=1,nDim-1 00403 m=AuxParFields(i) 00404 ip(m,2)=ip(m,2)+1 00405 iq(m,2)=iq(m,2)+1 00406 Info%aux(iq(1,1):iq(1,2),iq(2,1):iq(2,2),iq(3,1):iq(3,2),m) = & 00407 aux_psign*Info%aux(ip(1,1):ip(1,2),ip(2,1):ip(2,2),ip(3,1):ip(3,2),m) 00408 00409 ip(m,2)=ip(m,2)-1 00410 iq(m,2)=iq(m,2)-1 00411 00412 END DO 00413 iq(dim,:)=iq(dim,:)+(edge-1) 00414 ip(dim,:)=ip(dim,:)+(2-edge) 00415 Info%aux(iq(1,1):iq(1,2),iq(2,1):iq(2,2),iq(3,1):iq(3,2),dim) = & 00416 aux_nsign*Info%aux(ip(1,1):ip(1,2),ip(2,1):ip(2,2),ip(3,1):ip(3,2),dim) 00417 00418 END IF 00419 END DO 00420 CASE DEFAULT 00421 PRINT*, 'boundary condition type not implemented yet' 00422 STOP 00423 END SELECT 00424 END SUBROUTINE SetGhost 00425 END SUBROUTINE ApplyPhysicalBC 00426 00427 00430 SUBROUTINE CheckDivergence(Info, lStopDivergence) 00431 USE ProcessingDeclarations 00432 TYPE(InfoDef) :: Info 00433 INTEGER :: i 00434 REAL :: ri,rl,rh 00435 REAL(KIND=qPREC), DIMENSION(:,:,:), POINTER :: divergence 00436 LOGICAL, OPTIONAL :: lStopDivergence 00437 00438 IF (MaintainAuxArrays) THEN 00439 ALLOCATE(divergence(Info%mX(1),Info%mX(2),Info%mX(3))) 00440 00441 IF (nDim == 2) THEN 00442 IF (iCylindrical==NoCyl) THEN 00443 divergence(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3))=& 00444 (Info%aux(2:Info%mX(1)+1,1:Info%mX(2),1:Info%mX(3),1)+& 00445 Info%aux(1:Info%mX(1),2:Info%mX(2)+1,1:Info%mX(3),2)-& 00446 Info%aux(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),1) - & 00447 Info%aux(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),2))/levels(Info%level)%dx 00448 ELSE 00449 DO i=1,Info%mX(1) 00450 ri=1.d0/(Info%xBounds(1,1)+(REAL(i)-half)*levels(Info%level)%dx) 00451 rl=(Info%xBounds(1,1)+(i-1)*levels(Info%level)%dx) 00452 rh=(Info%xBounds(1,1)+(i)*levels(Info%level)%dx) 00453 divergence(i,1:Info%mX(2),1:Info%mX(3))=& 00454 (rl*ri*Info%aux(i+1,1:Info%mX(2),1:Info%mX(3),1)+& 00455 Info%aux(i,2:Info%mX(2)+1,1:Info%mX(3),2)-& 00456 rh*ri*Info%aux(i,1:Info%mX(2),1:Info%mX(3),1) - & 00457 Info%aux(i,1:Info%mX(2),1:Info%mX(3),2))/levels(Info%level)%dx 00458 END DO 00459 END IF 00460 ELSE 00461 divergence(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3))=& 00462 (Info%aux(2:Info%mX(1)+1,1:Info%mX(2),1:Info%mX(3),1)+& 00463 Info%aux(1:Info%mX(1),2:Info%mX(2)+1,1:Info%mX(3),2)+& 00464 Info%aux(1:Info%mX(1),1:Info%mX(2),2:Info%mX(3)+1,3)-& 00465 Info%aux(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),1) - & 00466 Info%aux(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),2) - & 00467 Info%aux(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),3))/levels(Info%level)%dx 00468 END IF 00469 00470 IF (ANY(ABS(divergence(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3))) > & 00471 1d-13/levels(Info%level)%dx*max(1d0, sqrt(SUM(Info%q(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),AuxFields(1:nDim))**2, 4))))) THEN 00472 ! CALL WriteDataFrame(100) 00473 write(*,*) "Divergence Warning", Info%level, MAXVAL(ABS(divergence(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3)))) 00474 00475 IF (PRESENT(lStopDivergence) .AND. lStopDivergence) THEN 00476 WRITE(*,*) "*** Simulation stopped due to divergence in the initialization! ***" 00477 STOP 00478 END IF 00479 END IF 00480 00481 DEALLOCATE(divergence) 00482 NULLIFY(divergence) 00483 END IF 00484 END SUBROUTINE CheckDivergence 00485 00490 SUBROUTINE JeansCriterionTest(rho, cs, errorflags, dx) 00491 00492 REAL(KIND=qPrec), DIMENSION(:,:,:) :: cs, rho 00493 INTEGER, DIMENSION(:,:,:) :: errorflags 00494 REAL(KIND=qPrec) :: jl_coeff, dx 00495 REAL(KIND=qPrec) :: jeans_length 00496 00497 jeans_length = dx * CellsPerJeansLength 00498 jl_coeff = PI / ScaleGrav 00499 00500 WHERE(SQRT(jl_coeff / rho) * cs < jeans_length) ErrorFlags = 1 00501 00502 ! rho > jl_coeff/(jeans_length)**2*(cs**2) 00503 00504 00505 ! JeansFact=sqrt(pi*gamma/ScaleGrav)/(JEAN_CELLS*levels(MaxLevel)%dx) 00506 !rhoJeans=JeansFact*sqrt(press(q(i,j,k,:))) Sink Particle test 00507 END SUBROUTINE JeansCriterionTest 00508 00509 END MODULE ModuleControl 00510