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 Christina_Original 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 !######################################################################### 00023 !============================================================================ 00024 ! BEARCLAW Astrophysical Applications 00025 !============================================================================ 00026 ! Christina Haig, Fabian Heitsch 00027 ! Department of Physics and Astronomy 00028 ! University of North Carolina, Chapel Hill 00029 ! Chapel Hill, NC 00030 !============================================================================ 00031 MODULE Problem 00032 ! User definitions of data structures associated with each node 00033 USE DataDeclarations 00034 USE CoolingSrc 00035 USE PhysicsDeclarations 00036 ! Global declarations for the BEARCLAW Astrophysical Applications 00037 USE GlobalDeclarations 00038 USE EOS 00039 USE Winds 00040 USE Ambients 00041 IMPLICIT NONE ! It's safer to require explicit declarations 00042 SAVE ! Save module information 00043 PRIVATE ! Everything is implicitly private, i.e. accessible only 00044 ! to this module. 00045 00046 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00047 ! ! 00048 ! Declarations Local to the Module ! 00049 ! ! 00050 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00051 00052 PUBLIC ProblemModuleInit, ProblemGridInit, CloudsBC, & 00053 ProblemBeforeStep, ProblemAfterStep, ProblemSetErrFlag, ProblemBeforeGlobalStep 00054 TYPE(CoolingDef),POINTER :: coolingobj 00055 INTEGER :: iCooling 00056 INTEGER :: iTracer1=0 00057 INTEGER :: iTracer2=0 00058 !PUBLIC :: scaleClouds,InitializeClouds,qinitClouds,b4stepClouds,afterstepClouds,CloudsBC!,CloudsErrFlag 00059 00060 !Cloud parameters 00061 REAL(KIND=qPrec) :: power,nCloud,TCloud,vx0,vy0,vz0,kx0,ky0,kz0,bmin,amaj,wght,lamda,sig0,pert,radius,slope 00062 REAL(KIND=xPrec) :: xlen,ylen,zlen,xmin,ymin,zmin 00063 INTEGER :: iseed, kmin, kmax 00064 REAL(KIND = qPrec),ALLOCATABLE,DIMENSION(:,:,:) :: amplitrand 00065 INTEGER :: mxtot,mytot,mztot 00066 00067 REAL(KIND=qPrec) :: alpha,beta ! analytic cooling parameters 00068 00069 NAMELIST/ProblemData/iCooling,iseed,power,kmin,kmax,nCloud,TCloud,vx0,vy0,vz0,kx0,ky0,kz0,bmin,amaj,sig0,pert,radius,slope 00070 00071 00072 CONTAINS 00073 00074 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00075 ! ! 00076 ! Public Module Prodecures ! 00077 ! ! 00078 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00079 00080 00081 SUBROUTINE ProblemModuleInit 00082 !REAL(KIND = qPrec),ALLOCATABLE,DIMENSION(:,:,:) :: amplitrand 00083 REAL(KIND=xPrec) :: dX(3), rhoOut, pOut, vxOut, vyOut, vzOut, BxOut, ByOut, BzOut 00084 REAL(KIND = qPrec) :: mina, maxa 00085 INTEGER :: iErr 00086 TYPE(pWindDef), DIMENSION(:), ALLOCATABLE :: MyWinds 00087 TYPE(AmbientDef), POINTER :: Ambient 00088 NAMELIST /AmbientData/ rhoOut, pOut, vxOut, vyOut, vzOut, BxOut, ByOut, BzOut 00089 00090 OPEN(UNIT=PROBLEM_DATA_HANDLE,IOSTAT=iErr,FILE='problem.data') 00091 READ(PROBLEM_DATA_HANDLE,NML=ProblemData) 00092 CALL CreateAmbient(Ambient) 00093 READ(PROBLEM_DATA_HANDLE,NML=AmbientData) 00094 Ambient%density=rhoOut 00095 Ambient%pressure=pOut 00096 Ambient%B(:)=(/BxOut, ByOut, BzOut/) 00097 Ambient%velocity(:)=(/vxOut, vyOut, vzOut/) 00098 00099 CLOSE(PROBLEM_DATA_HANDLE) 00100 00101 CALL AddTracer(iTracer1, 'FlowTracer_1') 00102 CALL AddTracer(iTracer2, 'FlowTracer_2') 00103 00104 mxtot = levels(MAXLEVEL)%mX(1);mytot=levels(MaxLEVEL)%mX(2);mztot=LEVELS(MaxLEVEL)%mX(3) 00105 xlen = GxBounds(1,2)-GxBounds(1,1) 00106 ylen = GxBounds(2,2)-GxBounds(2,1) 00107 zlen = GxBounds(3,2)-GxBounds(3,1) 00108 xmin = GxBounds(1,1) 00109 00110 ! This part sets the random interface perturbations 00111 IF (iSeed.NE.0) THEN 00112 00113 ALLOCATE(amplitrand(1:mytot,1:mztot,1)) 00114 CALL fldgen(mytot, mytot, 1, mztot, & 00115 mztot, mztot, 1, mztot, & 00116 1 , 1 , 1, 1 , & 00117 power , kmin , kmax , & 00118 iseed , amplitrand) 00119 00120 00121 mina = minval(amplitrand(1:mytot,1:mztot,1)) 00122 maxa = maxval(amplitrand(1:mytot,1:mztot,1)) 00123 amplitrand(:,:,:) = (pert*3.08568025d18/lscale)*((amplitrand(:,:,:)-mina)/(maxa-mina)-0.5d0) 00124 mina = minval(amplitrand(1:mytot,1:mztot,1)) 00125 maxa = maxval(amplitrand(1:mytot,1:mztot,1)) 00126 00127 END IF 00128 00129 00130 ! [BDS][20110106]: Do not create a new cooling object on a restart. The cooling object will be 00131 ! read in from the Chombo files. 00132 IF(iCooling>0) THEN 00133 IF (.NOT. lRestart) THEN 00134 ! see sources/cooling.f90::CreateCoolingObject for 00135 ! default values of a cooling source term 00136 CALL CreateCoolingObject(coolingobj) 00137 write(*,*) 'created cooling object' 00138 ELSE 00139 coolingobj => firstcoolingobj 00140 END IF 00141 END IF 00142 00143 coolingobj%iCooling=iCooling 00144 SELECT CASE(iCooling) ! cases defined in sources/cooling.f90 00145 CASE(NoCool) 00146 CASE(AnalyticCool) 00147 coolingobj%alpha=alpha 00148 coolingobj%beta=beta 00149 CASE(DMCool) 00150 CASE(IICool) 00151 CASE DEFAULT 00152 END SELECT 00153 00154 coolingobj%floortemp=0.001d0 00155 coolingobj%mintemp=0.001 00156 00157 00158 END SUBROUTINE ProblemModuleInit 00159 00160 SUBROUTINE ProblemGridInit(Info) 00161 00162 TYPE(InfoDef) :: Info ! Data associated with this grid 00163 00164 REAL (KIND = qPrec), POINTER, DIMENSION (:,:,:,:) :: q 00165 00166 ! The q-variable fields are as follows: 00167 ! 00168 ! q(:,:,:,:,1) = density 00169 ! q(:,:,:,:,2) = px 00170 ! q(:,:,:,:,3) = py 00171 ! q(:,:,:,:,4) = E 00172 ! q(:,:,:,:,5) = phi 00173 00174 INTEGER :: i, j, k, mx, my, mz,ix,iy,iz 00175 INTEGER :: mbc, rmbc, zrmbc 00176 REAL(KIND = xPrec) :: x, y, z, xl, yl, zl, dx, dy, dz, r,grdx,grdy,grdz,dx_maxlevel 00177 REAL(KIND = qPrec) :: mina, maxa, cosa, arlo, arhi,r2 00178 LOGICAL, SAVE :: lfirst = .TRUE. 00179 00180 00181 q => Info%q 00182 mbc = 0!Info%mbc !number of ghost cells 00183 rmbc=0 00184 00185 ! IF (Info%ErrorEstimation) THEN 00186 ! rmbc = (Info%AMRSteps-Info%AMRStep+1)*mbc 00187 ! ELSE 00188 ! rmbc = Info%r*mbc 00189 ! END IF 00190 !Erase this if not needed: mxtot = levels(MAXLEVEL)%mX(1) 00191 00192 mx = Info%mX(1);my = Info%mX(2);mz = Info%mX(3) 00193 dx = levels(Info%level)%dx;dy=dx;dz=dx !Info%dX(1);dy = Info%dX(2);dz = Info%dX(3) 00194 dx_maxlevel = levels(maxlevel)%dx 00195 xl = Info%XBounds(1,1);yl = Info%XBounds(2,1);zl = Info%XBounds(3,1) 00196 xmin = GxBounds(1,1);ymin = GxBounds(2,1);zmin = GxBounds(3,1) 00197 00198 q(:,:,:,:) = 0.d0 00199 00200 cosa = 1.d0 00201 00202 00203 SELECT CASE(nDim) 00204 CASE(2) 00205 zrmbc = 0; mz = 1; zl = 0; dz = 0 00206 00207 DO i = 1-rmbc, mx + rmbc 00208 x = (xl+(REAL(i,xPrec)-half)*dx) ! set x, y, z coordinates at cell centers 00209 grdx = (xl+REAL(i,xPrec)*dx) ! set x, y, z coordinates at cell walls to the right of cell centers 00210 ix = Info%mGlobal(1,1)+i-1 00211 DO j = 1, my 00212 y = (yl+(REAL(j,xPrec)-half)*dy) 00213 grdy = (yl+REAL(j,xPrec)*dy) 00214 iy = Info%mGlobal(2,1)+j-1 00215 00216 IF (iseed.EQ.0) THEN ! geometrical collision interface 00217 arlo = (xmin+ slope*y + 0.5d0*xlen - (grdx-dx))*dy & ! subtract everything below cell 00218 + (ylen*xlen*dabs(pert)/(2d0*pi*ky0)) & 00219 * (dsin(2d0*pi*ky0*grdy/ylen)-dsin(2d0*pi*ky0*(grdy-dy)/ylen)) 00220 00221 IF (kz0 .GT. 0d0) THEN ! this is for overlaying a perturbation at different scale. 00222 ! Scale is given by ky0/kz0. For ky0>kz0, the additional perturbation will have 00223 ! a larger scale at higher amplitude. Only for code stability tests. 00224 arlo = arlo + (ylen*xlen*dabs(pert*(ky0/kz0)**(5d0/3d0)) & 00225 /(2d0*pi*kz0))*(dsin(2d0*pi*kz0*grdy/ylen)-dsin(2d0*pi*kz0*(grdy-dy)/ylen)) 00226 END IF 00227 00228 arlo = arlo / (dx*dy) 00229 arhi = 1d0 - arlo 00230 ! Now we have the "lower" and "higher" area in each cell, i.e. the area (=volume...) 00231 ! of the current cell located left of perturbation (arlo) and right of perturbation (arhi). 00232 00233 ELSE ! Random interface 00234 00235 arlo = amplitrand(iy,1,1) - (grdx-dx) + slope*y 00236 arhi = grdx - amplitrand(iy,1,1) - slope*y 00237 00238 END IF 00239 00240 IF (arlo .LT. 0d0) THEN ! perturbation interface does not intersect cell. Cell to the right. 00241 arlo = 0d0 00242 arhi = 1d0 00243 END IF 00244 00245 IF (arhi .LT. 0d0) THEN ! perturbation interface does not intersect cell. Cell to the left. 00246 arhi = 0d0 00247 arlo = 1d0 00248 END IF 00249 00250 q(i,j,1,1) = nCloud 00251 q(i,j,1,2) = q(i,j,1,1)*(vx0*arlo - vx0*arhi)*cosa 00252 q(i,j,1,3) = q(i,j,1,1)*(vy0*arlo - vy0*arhi) 00253 00254 IF (radius.NE.0.d0) THEN ! we use a limited inflow region 00255 IF (dabs(y).LE.radius) THEN 00256 wght=1.d0 00257 ELSE 00258 wght=0.d0 00259 END IF 00260 ELSE 00261 wght=1.d0 00262 00263 END IF 00264 00265 q(i,j,1,2) = q(i,j,1,2) * wght 00266 00267 00268 IF (iTracer1 .GT.0) THEN 00269 IF (q(i,j,1,2).GT.1e-2) THEN 00270 q(i,j,1,iTracer1) = 1.d0 00271 ELSE 00272 q(i,j,1,iTracer1) = 1.d-10 00273 ENDIF 00274 END IF 00275 00276 IF (iTracer2.GT.0) THEN 00277 IF (q(i,j,1,2).LT.1e-2) THEN 00278 q(i,j,1,iTracer2) = 1.d0 00279 ELSE 00280 q(i,j,1,iTracer2) = 1.d-10 00281 ENDIF 00282 END IF 00283 00284 CALL CalcComp(q(i,j,1,:),TCloud) 00285 END DO 00286 END DO 00287 00288 00289 CASE(3) 00290 zrmbc = rmbc 00291 00292 DO i = 1-rmbc, mx + rmbc 00293 x = (xl+(REAL(i,xPrec)-half)*dx) ! set x, y, z coordinates at cell centers 00294 grdx = (xl+REAL(i,xPrec)*dx) ! set x, y, z coordinates at cell walls to the right of cell centers 00295 !ix = Info%mGlobal(1,1)+i-1 00296 ix = NInt(half+(x-xmin)/dx_maxlevel) 00297 DO j = 1, my 00298 y = (yl+(REAL(j,xPrec)-half)*dy) 00299 grdy = (yl+REAL(j,xPrec)*dy) 00300 !iy = Info%mGlobal(2,1)+j-1 00301 iy = NInt(half+(y-ymin)/dx_maxlevel) 00302 DO k = 1, mz 00303 z = (zl+(REAL(k,xPrec)-half)*dz) 00304 grdz = (zl+REAL(k,xPrec)*dz) 00305 !iz = Info%mGlobal(3,1)+k-1 00306 iz = NInt(half+(z-zmin/dx_maxlevel)) 00307 00308 00309 IF (iseed.EQ.0) THEN ! geometrical collision interface 00310 arlo = (xmin+0.5d0*xlen-(grdx-dx)+slope*y)*dy*dz + (zlen*ylen*xlen*dabs(pert)/(4d0*pi**2*ky0*kz0))*& 00311 (dsin(2d0*pi*ky0*grdy/ylen)-dsin(2d0*pi*ky0*(grdy-dy)/ylen)) & 00312 * (dsin(2d0*pi*kz0*grdz/zlen)-dsin(2d0*pi*kz0*(grdz-dz)/zlen)) 00313 arlo = arlo / (dx*dy*dz) 00314 arhi = 1d0 - arlo 00315 ELSE ! random collision interface. Needs to be centered on x = 0d0. 00316 00317 00318 arlo = amplitrand(iy,iz,1) - (grdx-dx) + slope*y 00319 arhi = grdx - amplitrand(iy,iz,1) - slope*y 00320 00321 END IF 00322 00323 !IF (radius .GT. 0d0) THEN ! mimicking taper due to expanding blast wave 00324 ! cosa = dcos(datan(dsqrt(y**2 + z**2)/radius)) 00325 !END IF 00326 00327 IF (arlo .LT. 0d0) THEN ! perturbation interface does not intersect cell. Cell to the right. 00328 arlo = 0d0 00329 arhi = 1d0 00330 END IF 00331 00332 IF (arhi .LT. 0d0) THEN ! perturbation interface does not intersect cell. Cell to the left. 00333 arhi = 0d0 00334 arlo = 1d0 00335 END IF 00336 00337 q(i,j,k,1) = nCloud 00338 q(i,j,k,2) = nCloud*(vx0*arlo - vx0*arhi)*cosa 00339 q(i,j,k,3) = nCloud*(vy0*arlo - vy0*arhi) 00340 q(i,j,k,4) = nCloud*(vz0*arlo - vz0*arhi) 00341 00342 wght=1.d0 00343 00344 IF ((bmin .NE. 0d0) .AND. (amaj .NE. 0d0)) THEN ! we use a limited inflow region 00345 r2 = ((y/amaj)**2+(z/bmin)**2) ! normalized elliptical "radius" 00346 00347 IF (sig0 .EQ. 0d0) THEN ! box profile 00348 IF (r2 .GT. 1d0) THEN ! ellipse condition 00349 wght = 0d0 ! outside inflow region 00350 ELSE 00351 wght = 1d0 ! inside inflow region 00352 END IF 00353 ELSE IF (sig0 .LT. 0d0) THEN ! Gaussian profile 00354 wght = dexp(-r2/(sig0/xlen)**2) ! needs to be normalized 00355 ELSE IF (sig0 .GT. 0d0) THEN ! double tanh profile as a taper for viscous profile (BGK!) 00356 r2 = dsqrt(r2) 00357 wght = 0.5d0*(1d0-(dexp((r2-1d0)*xlen/sig0)-dexp(-(r2-1d0)*xlen/sig0))& 00358 /(dexp((r2-1d0)*xlen/sig0)+dexp(-(r2-1d0)*xlen/sig0))) 00359 END IF 00360 00361 q(i,j,k,2) = q(i,j,k,2) * wght 00362 END IF 00363 00364 00365 IF (iTracer1.GT.0) THEN 00366 IF ((q(i,j,k,2).GT.0).AND.(wght.GT.1e-2)) THEN 00367 q(i,j,k,iTracer1) = 1.d0 00368 ELSE 00369 q(i,j,k,iTracer1) = 1.d-10 00370 ENDIF 00371 END IF 00372 00373 IF (iTracer2.GT.0) THEN 00374 IF ((q(i,j,k,2).LT.0).AND.(wght.GT.1e-2)) THEN 00375 q(i,j,k,iTracer2) = 1.d0 00376 ELSE 00377 q(i,j,k,iTracer2) = 1.d-10 00378 ENDIF 00379 END IF 00380 00381 CALL CalcComp(q(i,j,k,:),TCloud) 00382 00383 00384 END DO 00385 END DO 00386 END DO 00387 00388 CASE DEFAULT 00389 PRINT *, "Error in qinitCloud, invalid problem dimensions" 00390 STOP 00391 END SELECT 00392 00393 00394 00395 END SUBROUTINE ProblemGridInit 00396 00397 SUBROUTINE ProblemBeforeSTep(Info) 00398 TYPE(InfoDef) :: Info 00399 CALL CloudsBC(Info) 00400 END SUBROUTINE ProblemBeforeSTep 00401 00402 SUBROUTINE ProblemSetErrFlag(Info) 00403 TYPE(InfoDef) :: Info 00404 END SUBROUTINE ProblemSetErrFlag 00405 00406 SUBROUTINE ProblemAfterStep(Info) 00407 TYPE(InfoDef) :: Info 00408 REAL (KIND = qPrec), POINTER, DIMENSION (:,:,:,:) :: q 00409 00410 INTEGER :: i, j, k, mx, my, mz 00411 INTEGER :: mbc, rmbc, zrmbc, ibc ! See assignment statement 00412 REAL(KIND = xPrec) :: x, y, z, xl, yl, zl, dx, dy, dz 00413 REAL(KIND = qPrec) :: Temp 00414 REAL (KIND=qPrec), PARAMETER :: F = 2d-26 00415 LOGICAL, SAVE :: lfirst = .TRUE. 00416 00417 q => Info%q 00418 mbc = 0!Info%mbc ! number of ghost cells 00419 rmbc=0 00420 00421 ! IF (Info%ErrorEstimation) THEN 00422 ! rmbc = (Info%AMRSteps-Info%AMRStep+1)*mbc 00423 ! ELSE 00424 ! rmbc = Info%r*mbc 00425 ! END IF 00426 00427 00428 mx = Info%mX(1);my = Info%mX(2);mz = Info%mX(3) 00429 dx = levels(Info%level)%dx;dy=dx;dz=dx!;Info%dX(1);dy = Info%dX(2);dz = Info%dX(3) 00430 xl = Info%XBounds(1,1);yl = Info%XBounds(2,1);zl = Info%xBounds(3,1) 00431 00432 SELECT CASE(nDim) 00433 CASE(2) 00434 zrmbc = 0; mz = 1; zl = 0; dz = 0 00435 00436 DO i = 1-rmbc, mx + rmbc 00437 x = (xl+(REAL(i,xPrec)-half)*dx) ! set x, y coordinates at cell centers 00438 DO j = 1-rmbc, my + rmbc 00439 y = (yl+(REAL(j,xPrec)-half)*dy) 00440 00441 ! Floor on both metallicity tracers 00442 IF (iTracer1.GT.0) THEN 00443 IF (q(i,j,1,iTracer1).LT.1.d-10) q(i,j,1,iTracer1)=1.d-10 00444 END IF 00445 00446 IF (iTracer2.GT.0) THEN 00447 IF (q(i,j,1,iTracer2).LT.1.d-10) q(i,j,1,iTracer2)=1.d-10 00448 END IF 00449 00450 00451 ! Temp = Tempscale*Press(q(i,j,1,:))/q(i,j,1,1) 00452 ! q(i,j,1,iTracer1) = q(i,j,1,1)*nscale*F*(1-q(i,j,1,1)*nscale*IICoolingRate(Temp)) 00453 !IF (i==1.AND.j==1) PRINT*,IIColingRate(Temp) 00454 END DO 00455 END DO 00456 00457 00458 00459 00460 CASE(3) 00461 zrmbc = rmbc 00462 DO i = 1-rmbc, mx + rmbc 00463 x = (xl+(REAL(i,xPrec)-half)*dx) ! set x, y, z coordinates at cell centers 00464 DO j = 1-rmbc, my + rmbc 00465 y = (yl+(REAL(j,xPrec)-half)*dy) 00466 DO k = 1-rmbc, mz + rmbc 00467 z = (zl+(REAL(k,xPrec)-half)*dz) 00468 00469 ! Floor on both metallicity tracers 00470 IF (iTracer1.GT.0) THEN 00471 IF (q(i,j,k,iTracer1).LT.1.d-10) q(i,j,k,iTracer1)=1.d-10 00472 END IF 00473 00474 IF (iTracer2.GT.0) THEN 00475 IF (q(i,j,k,iTracer2).LT.1.d-10) q(i,j,k,iTracer2)=1.d-10 00476 END IF 00477 00478 00479 ! Temp = Tempscale*Press(q(i,j,1,:))/q(i,j,1,1) 00480 ! q(i,j,1,iTracer1) = q(i,j,1,1)*nscale*F*(1-q(i,j,1,1)*nscale*IICoolingRate(Temp)) 00481 !IF (i==1.AND.j==1) PRINT*,IICoolingRate(Temp) 00482 END DO 00483 END DO 00484 END DO 00485 00486 00487 CASE DEFAULT 00488 PRINT *, "Error in b4stepCloud, invalid problem dimensions" 00489 STOP 00490 END SELECT 00491 00492 00493 END SUBROUTINE ProblemAfterStep 00494 00495 00496 00497 SUBROUTINE CloudsBC(Info) 00498 00499 USE EOS 00500 TYPE(InfoDef) :: Info 00501 00502 REAL (KIND = qPrec), POINTER, DIMENSION (:,:,:,:) :: q 00503 00504 INTEGER :: i, j, k, mx, my, mz,ibc 00505 INTEGER :: mbc, rmbc, zrmbc ! See assignment statement 00506 REAL(KIND = xPrec) :: x, y, z, xl, yl, zl, dx, dy, dz, r2 00507 REAL(KIND = qPrec) :: Temp 00508 REAL (KIND=qPrec), PARAMETER :: F = 2d-26 00509 LOGICAL, SAVE :: lfirst = .TRUE. 00510 00511 q => Info%q 00512 mbc = levels(Info%level)%ambc(levels(Info%level)%step) ! number of ghost cells 00513 rmbc = mbc 00514 mx = Info%mX(1);my = Info%mX(2);mz = Info%mX(3) 00515 dx = levels(Info%level)%dx !Info%dX(1);dy = Info%dX(2);dz = Info%dX(3) 00516 dy=dx 00517 dz=dx 00518 xl = Info%XBounds(1,1);yl = Info%XBounds(2,1);zl = Info%XBounds(3,1) 00519 00520 00521 SELECT CASE(nDim) 00522 CASE(2) 00523 zrmbc = 0; mz = 1; zl = 0; dz = 0 00524 00525 IF (Info%xBounds(1,1)==GxBounds(1,1) .AND. Gmthbc(1,1)==10) THEN ! Left Side 00526 00527 DO i = 1-rmbc, 0 00528 !x = (xl+(REAL(i,xPrec)-half)*dx) ! set x, y coordinates at cell centers 00529 DO j = 1-rmbc, my+rmbc 00530 y = (yl+(REAL(j,xPrec)-half)*dy) 00531 00532 ! First set extrapolating conditions, just in case 00533 Info%q(i,j,1,:)=Info%q(1,j,1,:) 00534 00535 IF (radius.NE.0.d0) THEN ! we use a limited inflow region 00536 00537 IF (dabs(y).LE.radius) THEN ! Inside 'circle' 00538 q(i,j,1,1) = nCloud 00539 q(i,j,1,2) = q(i,j,1,1)*vx0 00540 q(i,j,1,3) = q(i,j,1,1)*vy0 00541 CALL CalcComp(q(i,j,1,:),TCloud) 00542 00543 IF (iTracer1.GT.0) THEN 00544 q(i,j,1,iTracer1) = 1.d0 00545 END IF 00546 00547 IF (iTracer2.GT.0) THEN 00548 q(i,j,1,iTracer2) = 1.d-10 00549 END IF 00550 00551 ELSE ! Outside 'circle' 00552 IF (q(i,j,1,2).GT.0.d0) THEN 00553 Temp = TempScale*Press(q(i,j,1,:))/q(i,j,1,1) 00554 q(i,j,1,2)=0.d0 00555 CALL CalcComp(q(i,j,1,:),Temp) 00556 END IF 00557 00558 END IF 00559 00560 ELSE ! Inflow region is not limited (entire left side) 00561 q(i,j,1,1) = nCloud 00562 q(i,j,1,2) = q(i,j,1,1)*vx0 00563 q(i,j,1,3) = q(i,j,1,1)*vy0 00564 CALL CalcComp(q(i,j,1,:),TCloud) 00565 00566 IF (iTracer1.GT.0) THEN 00567 q(i,j,1,iTracer1) = 1.d0 00568 END IF 00569 00570 IF (iTracer2.GT.0) THEN 00571 q(i,j,1,iTracer2) = 1.d-10 00572 END IF 00573 00574 END IF 00575 00576 END DO 00577 END DO 00578 00579 END IF 00580 00581 IF (Info%xBounds(1,2)==GxBounds(1,2) .AND. Gmthbc(1,2)==10) THEN ! Right Side 00582 00583 DO i = mx+1, mx+rmbc 00584 !x = (xl+(REAL(i,xPrec)-half)*dx) ! set x, y coordinates at cell centers 00585 DO j = 1-rmbc, my+rmbc 00586 y = (yl+(REAL(j,xPrec)-half)*dy) 00587 00588 ! First set boundary conditions to extrapolating, just in case 00589 Info%q(i,j,1,:)=Info%q(mx,j,1,:) 00590 00591 IF (radius.NE.0.d0) THEN ! we use a limited inflow region 00592 00593 IF (dabs(y).LE.radius) THEN ! Inside 'circle' 00594 q(i,j,1,1) = nCloud 00595 q(i,j,1,2) = -q(i,j,1,1)*vx0 00596 q(i,j,1,3) = q(i,j,1,1)*vy0 00597 CALL CalcComp(q(i,j,1,:),TCloud) 00598 00599 IF (iTracer1.GT.0) THEN 00600 q(i,j,1,iTracer1) = 1.d-10 00601 END IF 00602 IF (iTracer2.GT.0) THEN 00603 q(i,j,1,iTracer2) = 1.d0 00604 END IF 00605 00606 ELSE ! Outside 'circle' 00607 IF (q(i,j,1,2).LT.0.d0) THEN 00608 Temp = TempScale*Press(q(i,j,1,:))/q(i,j,1,1) 00609 q(i,j,1,2)=0.d0 00610 CALL CalcComp(q(i,j,1,:),Temp) 00611 END IF 00612 END IF 00613 00614 ELSE ! Inflow region is not limited (entire right side) 00615 q(i,j,1,1) = nCloud 00616 q(i,j,1,2) = -q(i,j,1,1)*vx0 00617 q(i,j,1,3) = q(i,j,1,1)*vy0 00618 CALL CalcComp(q(i,j,1,:),TCloud) 00619 00620 IF (iTracer1.GT.0) THEN 00621 q(i,j,1,iTracer1) = 1.d-10 00622 END IF 00623 IF (iTracer2.GT.0) THEN 00624 q(i,j,1,iTracer2) = 1.d0 00625 END IF 00626 00627 END IF 00628 00629 00630 END DO 00631 END DO 00632 END IF 00633 00634 00635 IF (Info%xBounds(2,1)==GxBounds(2,1) .AND. Gmthbc(2,1)==10) THEN ! Bottom 00636 00637 DO ibc=1,rmbc 00638 Info%q(:,1-ibc,1,:) = Info%q(:,1,1,:) 00639 END DO 00640 00641 DO i = 1-rmbc, mx+rmbc 00642 !x = (xl+(REAL(i,xPrec)-half)*dx) ! set x, y coordinates at cell centers 00643 DO j = 1-rmbc, 0 00644 !y = (yl+(REAL(j,xPrec)-half)*dy) 00645 00646 ! Extrapolating 00647 Info%q(i,j,1,:) = Info%q(i,1,1,:) 00648 00649 IF (q(i,j,1,3).GT.0.d0) THEN 00650 Temp = TempScale*Press(q(i,j,1,:))/q(i,j,1,1) 00651 q(i,j,1,3)=0.d0 00652 CALL CalcComp(q(i,j,1,:),Temp) 00653 END IF 00654 00655 00656 END DO 00657 END DO 00658 00659 END IF 00660 00661 IF (Info%xBounds(2,2)==GxBounds(2,2) .AND. Gmthbc(2,2)==10) THEN ! Top 00662 00663 00664 DO i = 1-rmbc, mx+rmbc 00665 !x = (xl+(REAL(i,xPrec)-half)*dx) ! set x, y coordinates at cell centers 00666 DO j = my+1, my+rmbc 00667 !y = (yl+(REAL(j,xPrec)-half)*dy) 00668 00669 ! Extrapolating 00670 Info%q(i,j,1,:) = Info%q(i,my,1,:) 00671 00672 IF (q(i,j,1,3).LT.0.d0) THEN 00673 00674 Temp = TempScale*Press(q(i,j,1,:))/q(i,j,1,1) 00675 q(i,j,1,3)=0.d0 00676 CALL CalcComp(q(i,j,1,:),Temp) 00677 END IF 00678 END DO 00679 END DO 00680 00681 END IF 00682 00683 CASE(3) 00684 zrmbc = rmbc 00685 00686 00687 IF (Info%xBounds(1,1)==GxBounds(1,1) .AND. Gmthbc(1,1)==10) THEN ! Left Side 00688 00689 DO i = 1-rmbc, 0 00690 !x = (xl+(REAL(i,xPrec)-half)*dx) ! set x, y coordinates at cell centers 00691 DO j = 1-rmbc, my+rmbc 00692 y = (yl+(REAL(j,xPrec)-half)*dy) 00693 DO k = 1-rmbc, mz+rmbc 00694 z = (zl+(REAL(k,xPrec)-half)*dz) 00695 00696 ! First set extrapolating conditions, just in case 00697 Info%q(i,j,k,:)=Info%q(1,j,k,:) 00698 00699 IF ((bmin .NE. 0d0) .AND. (amaj .NE. 0d0)) THEN ! we use a limited inflow region 00700 r2 = ((y/amaj)**2+(z/bmin)**2) ! normalized elliptical "radius" 00701 00702 IF (sig0 .EQ. 0d0) THEN ! box profile 00703 IF (r2 .GT. 1d0) THEN ! ellipse condition 00704 wght = 0d0 ! outside inflow region 00705 ELSE 00706 wght = 1d0 ! inside inflow region 00707 END IF 00708 00709 ELSE IF (sig0 .LT. 0d0) THEN ! Gaussian profile 00710 wght = dexp(-r2/(sig0/xlen)**2) ! needs to be normalized 00711 00712 ELSE IF (sig0 .GT. 0d0) THEN ! double tanh profile as a taper for viscous profile (BGK!) 00713 r2 = dsqrt(r2) 00714 wght = 0.5d0*(1d0-(dexp((r2-1d0)*xlen/sig0)-dexp(-(r2-1d0)*xlen/sig0))& 00715 /(dexp((r2-1d0)*xlen/sig0)+dexp(-(r2-1d0)*xlen/sig0))) 00716 END IF 00717 00718 IF (wght.gt.1e-2) THEN ! Consider yourself in inflow region 00719 q(i,j,k,1) = nCloud 00720 q(i,j,k,2) = q(i,j,k,1)*vx0*wght 00721 q(i,j,k,3) = q(i,j,k,1)*vy0*wght 00722 CALL CalcComp(q(i,j,k,:),TCloud) 00723 IF (iTracer1.GT.0) q(i,j,k,iTracer1) = 1.d0 00724 IF (iTracer2.GT.0) q(i,j,k,ITracer1+1) = 1.d-10 00725 00726 ELSE ! Outflow conditions only outside of inflow region 00727 IF (q(i,j,k,2).GT.0.d0) THEN 00728 Temp = TempScale*Press(q(i,j,k,:))/q(i,j,k,1) 00729 q(i,j,k,2)=0.d0 00730 CALL CalcComp(q(i,j,k,:),Temp) 00731 END IF 00732 00733 END IF 00734 00735 ELSE ! Inflow is not limited - entire left side 00736 00737 q(i,j,k,1) = nCloud 00738 q(i,j,k,2) = q(i,j,k,1)*vx0 00739 q(i,j,k,3) = q(i,j,k,1)*vy0 00740 CALL CalcComp(q(i,j,k,:),TCloud) 00741 00742 IF (iTracer1.GT.0) q(i,j,k,iTracer1) = 1.d0 00743 IF (iTracer2.GT.0) q(i,j,k,iTracer2) = 1.d-10 00744 00745 END IF 00746 00747 END DO 00748 END DO 00749 END DO 00750 00751 END IF 00752 IF (Info%xBounds(1,2)==GxBounds(1,2) .AND. Gmthbc(1,2)==10) THEN ! Right Side 00753 00754 DO i = mx+1, mx+rmbc 00755 !x = (xl+(REAL(i,xPrec)-half)*dx) ! set x, y coordinates at cell centers 00756 DO j = 1-rmbc, my+rmbc 00757 y = (yl+(REAL(j,xPrec)-half)*dy) 00758 DO k = 1-rmbc, mz+rmbc 00759 z = (zl+(REAL(k,xPrec)-half)*dz) 00760 00761 ! First set extrapolating conditions, just in case 00762 Info%q(i,j,k,:)=Info%q(mx,j,k,:) 00763 00764 IF ((bmin .NE. 0d0) .AND. (amaj .NE. 0d0)) THEN ! we use a limited inflow region 00765 r2 = ((y/amaj)**2+(z/bmin)**2) ! normalized elliptical "radius" 00766 00767 IF (sig0 .EQ. 0d0) THEN ! box profile 00768 IF (r2 .GT. 1d0) THEN ! ellipse condition 00769 wght = 0d0 ! outside inflow region 00770 ELSE 00771 wght = 1d0 ! inside inflow region 00772 END IF 00773 00774 ELSE IF (sig0 .LT. 0d0) THEN ! Gaussian profile 00775 wght = dexp(-r2/(sig0/xlen)**2) ! needs to be normalized 00776 00777 ELSE IF (sig0 .GT. 0d0) THEN ! double tanh profile as a taper for viscous profile (BGK!) 00778 r2 = dsqrt(r2) 00779 wght = 0.5d0*(1d0-(dexp((r2-1d0)*xlen/sig0)-dexp(-(r2-1d0)*xlen/sig0))& 00780 /(dexp((r2-1d0)*xlen/sig0)+dexp(-(r2-1d0)*xlen/sig0))) 00781 END IF 00782 00783 IF (wght.gt.1e-2) THEN 00784 q(i,j,k,1) = nCloud 00785 q(i,j,k,2) = -q(i,j,k,1)*vx0*wght 00786 q(i,j,k,3) = -q(i,j,k,1)*vy0*wght 00787 CALL CalcComp(q(i,j,k,:),TCloud) 00788 IF (iTracer1.GT.0) q(i,j,k,iTracer1) = 1.d-10 00789 IF (iTracer2.GT.0) q(i,j,k,ITracer1+1) = 1.d0 00790 00791 ELSE 00792 IF (q(i,j,k,2).LT.0.d0) THEN 00793 Temp = TempScale*Press(q(i,j,k,:))/q(i,j,k,1) 00794 q(i,j,k,2)=0.d0 00795 CALL CalcComp(q(i,j,k,:),Temp) 00796 END IF 00797 00798 END IF 00799 00800 ELSE ! Inflow is not limited - entire right side 00801 00802 q(i,j,k,1) = nCloud 00803 q(i,j,k,2) = -q(i,j,k,1)*vx0 00804 q(i,j,k,3) = q(i,j,k,1)*vy0 00805 CALL CalcComp(q(i,j,k,:),TCloud) 00806 00807 IF (iTracer1.GT.0) q(i,j,k,iTracer1) = 1.d-10 00808 IF (iTracer2.GT.0) q(i,j,k,iTracer2) = 1.d0 00809 00810 END IF 00811 00812 END DO 00813 END DO 00814 END DO 00815 00816 END IF 00817 IF (Info%xBounds(2,1)==GxBounds(2,1) .AND. Gmthbc(2,1)==10) THEN ! bottom 00818 00819 DO i = 1-rmbc, mx+rmbc 00820 !x = (xl+(REAL(i,xPrec)-half)*dx) ! set x, y coordinates at cell centers 00821 DO j = 1-rmbc, 0 00822 !y = (yl+(REAL(j,xPrec)-half)*dy) 00823 DO k = 1-rmbc, mz+rmbc 00824 00825 ! Extrapolating 00826 Info%q(i,j,k,:) = Info%q(i,k,1,:) 00827 00828 IF (q(i,j,k,3).GT.0.d0) THEN 00829 Temp = TempScale*Press(q(i,j,k,:))/q(i,j,k,1) 00830 q(i,j,k,3)=0.d0 00831 CALL CalcComp(q(i,j,k,:),Temp) 00832 END IF 00833 00834 00835 END DO 00836 END DO 00837 END DO 00838 00839 END IF 00840 00841 IF (Info%xBounds(2,2)==GxBounds(2,2) .AND. Gmthbc(2,2)==10) THEN ! Top Side 00842 00843 DO i = 1-rmbc, mx+rmbc 00844 !x = (xl+(REAL(i,xPrec)-half)*dx) ! set x, y coordinates at cell centers 00845 DO j = my+1, my+rmbc 00846 !y = (yl+(REAL(j,xPrec)-half)*dy) 00847 DO k = 1-rmbc, mz+rmbc 00848 00849 ! Extrapolating 00850 Info%q(i,j,k,:) = Info%q(i,my,k,:) 00851 00852 IF (q(i,j,k,3).LT.0.d0) THEN 00853 Temp = TempScale*Press(q(i,j,k,:))/q(i,j,k,1) 00854 q(i,j,k,3)=0.d0 00855 CALL CalcComp(q(i,j,k,:),Temp) 00856 END IF 00857 00858 END DO 00859 END DO 00860 END DO 00861 00862 END IF 00863 00864 IF (Info%xBounds(3,1)==GxBounds(3,1) .AND. Gmthbc(3,1)==10) THEN ! front Side 00865 00866 DO i = 1-rmbc, mx+rmbc 00867 !x = (xl+(REAL(i,xPrec)-half)*dx) ! set x, y coordinates at cell centers 00868 DO j = 1-rmbc, my+rmbc 00869 !y = (yl+(REAL(j,xPrec)-half)*dy) 00870 DO k = 1-rmbc, 0 00871 00872 ! Extrapolating 00873 Info%q(i,j,k,:) = Info%q(i,j,1,:) 00874 00875 IF (q(i,j,k,4).GT.0.d0) THEN 00876 Temp = TempScale*Press(q(i,j,k,:))/q(i,j,k,1) 00877 q(i,j,k,4)=0.d0 00878 CALL CalcComp(q(i,j,k,:),Temp) 00879 END IF 00880 00881 END DO 00882 END DO 00883 END DO 00884 00885 END IF 00886 00887 IF (Info%xBounds(3,2)==GxBounds(3,2) .AND. Gmthbc(3,2)==10) THEN ! back Side 00888 00889 00890 DO i = 1-rmbc, mx+rmbc 00891 !x = (xl+(REAL(i,xPrec)-half)*dx) ! set x, y coordinates at cell centers 00892 DO j = 1-rmbc, my+rmbc 00893 !y = (yl+(REAL(j,xPrec)-half)*dy) 00894 DO k = mz+1, mz+rmbc 00895 00896 ! Extrapolating 00897 Info%q(i,j,k,:) = Info%q(i,j,mz,:) 00898 00899 IF (q(i,j,k,4).LT.0.d0) THEN 00900 00901 Temp = TempScale*Press(q(i,j,k,:))/q(i,j,k,1) 00902 q(i,j,k,4)=0.d0 00903 CALL CalcComp(q(i,j,k,:),Temp) 00904 END IF 00905 END DO 00906 END DO 00907 END DO 00908 00909 END IF 00910 00911 00912 CASE DEFAULT 00913 PRINT *, "Error in CloudsBC, invalid problem dimensions" 00914 STOP 00915 END SELECT 00916 00917 00918 00919 00920 END SUBROUTINE CloudsBC 00921 00922 00923 SUBROUTINE CalcComp(q,T,P) 00924 REAL(KIND=xprec), DIMENSION(nRvars) :: q 00925 REAL(KIND=xprec),INTENT(IN) :: T 00926 REAL(KIND=xprec),INTENT(OUT),OPTIONAL :: P 00927 REAL(KIND=xprec) :: gam,mu,Tin,Pin 00928 ! REAL (KIND=qPrec), DIMENSION(0:nSpeciesHI) :: nvec 00929 00930 mu=Xmu 00931 gam=gamma 00932 00933 ! IF(icooling==2) THEN 00934 ! for chemical flow, get the molecular weight and 00935 ! gamma for computing energy 00936 ! CALL EOS_vars(q,nvec=nvec(iHaux:nSpeciesHI),mu=mu,gamma=gam) 00937 ! ELSE IF(iEOS==1) THEN 00938 ! gam=gammac 00939 ! END IF 00940 00941 ! T is scaled by TempScale from Physics.data, so input in real units 00942 q(iE)= q(1)*one/(gam-one)*T/(TempScale*mu) 00943 00944 Pin=(gam-one)*q(iE) 00945 IF(PRESENT(P)) P=Pin 00946 00947 q(iE) = q(iE)+half*DOT_PRODUCT(q(2:m_high),q(2:m_high))/q(1) 00948 00949 ! IF(iEntropy/=0) THEN 00950 ! q(iEntropy) = Pin*q(1)**(one-gam) 00951 ! END IF 00952 END SUBROUTINE CalcComp 00953 ! End CalcComp 00954 00955 INCLUDE 'fldgen.f90' 00956 INCLUDE 'genak.f90' 00957 INCLUDE 'nr_fourn.f90' 00958 INCLUDE 'nr_randnum.f90' 00959 00960 SUBROUTINE ProblemBeforeGlobalStep(n) 00961 INTEGER :: n 00962 END SUBROUTINE ProblemBeforeGlobalStep 00963 00964 00965 END MODULE Problem