Scrambler  1
GravitationalCascade/problem.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 !    problem.f90 of module GravitationalCascade 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 MODULE Problem
00036   USE DataDeclarations
00037   USE Clumps
00038   USE ProcessingDeclarations  
00039   USE CoolingSrc
00040   USE Winds  
00041   USE Ambients
00042   USE Histograms
00043   USE Fields
00044   USE Totals
00045   USE PDFs
00046   USE Projections
00047   USE Refinements
00048   IMPLICIT NONE
00049   SAVE
00050 
00051   PUBLIC ProblemModuleInit, ProblemGridInit, &
00052        ProblemBeforeStep, ProblemAfterStep, ProblemSetErrFlag, ProblemBeforeGlobalStep
00053   
00054   REAL(KIND=qPREC) :: r_children(0:MaxDepth), rho_children(0:MaxDepth), fill_fraction(0:MaxDepth), xi(0:MaxDepth)
00055   INTEGER :: ClumpTreeDepth=1, n_children
00056 !   TYPE ClumpTreeDef
00057 !      TYPE(Clumpdef), POINTER :: Clump
00058 !      TYPE(ClumpTreeListDef), POINTER :: Children      
00059 !   END type ClumpTreeDef!
00060 
00061 !   TYPE ClumpTreeListDef
00062 !      TYPE(ClumpTreeDef), POINTER :: self
00063 !      TYPE(ClumpTreeDef), POINTER :: next
00064 !   END type ClumpTreeListDef
00065 
00066 CONTAINS
00067 
00068 
00070    SUBROUTINE ProblemModuleInit()      
00071       TYPE(CoolingDef),POINTER :: coolingobj
00072       TYPE(AmbientDef), POINTER :: Ambient
00073       TYPE(ClumpDef), POINTER :: Clump
00074       TYPE(HistogramDef), POINTER :: HISTOGRAM
00075       TYPE(PDFDef), POINTER :: PDF
00076       TYPE(ProjectionDef), POINTER :: Projection
00077       TYPE(WindDef), POINTER :: Wind
00078       REAL(KIND=qPREC) :: A,B,ValidRoots(3),r_cloud,r_ambient,chi
00079       INTEGER :: i,j,nroots,  edge, dim, n
00080       COMPLEX(8) :: MyRoots(3)
00081       Logical :: lUseExisting = .False.     
00082 
00083       NAMELIST /ProblemData/ &
00084            r_cloud, &              !Radius of single cloud
00085            r_ambient, &            !Effective radius of ambient cloud
00086            n_children, &           !number of child clumps
00087            chi, &                  !desired density contrast between child and parent
00088            ClumpTreeDepth, &          !Depth of nested structures  !0 is no nesting
00089            lUseExisting            ! flag for checking if clumplets.data exist
00090       OPEN(UNIT=PROBLEM_DATA_HANDLE, FILE='problem.data', STATUS="OLD")
00091       READ(PROBLEM_DATA_HANDLE,NML=ProblemData)
00092 
00093       ! Begin with initial size and recursively use the number of child clumps and the filling fractions to generate children of sizes...  Then set the densities of the finest children and work backwards to set the densities of the parents to give the desired mean densities of the parent.  This does however set a limit on the filling fraction?  For example imagine 7 clumps of r=r_parent/3 for a filling fraction of 7/27?  Since the Jeans length scales like rho^-1/2, the "jeans density" scales like r^{-2} so densities of these child clumps would need to be 9 times that of the parent jeans density.  However for the parent cloud to have a mean density 1/9th of the child densities, the child filling fraction must be less then 1/9th so only 3 of these clouds could exist for the parent 'cloud' to be appropriately bound...  (albeit with a density of 0).
00094 
00095 
00096       xi(ClumpTreeDepth)=1d0
00097       DO n=ClumpTreeDepth-1,0, -1
00098          B=1d0/(xi(n+1)*chi)
00099          A=1d0-B
00100          MyRoots=GetRoots(n_children**2*A**3,3d0*n_children**2*A**2*B-1d0,3*n_children**2*A*B**2,n_children**2*B**3,nroots)
00101          j=0
00102          DO i=1,nroots 
00103             IF (AIMAG(MyRoots(i)) < tiny(1d0) .AND. REAL(MyRoots(i)) > tiny(1d0) .AND. REAL(1d0-MyRoots(i)) > tiny(1d0)) THEN
00104                j=j+1
00105                ValidRoots(j)=REAL(MyRoots(i))
00106             END IF
00107          END DO
00108          write(*,*) myroots(1:nroots)
00109          !         write(*,*) n_children**2*A**3*myroots**3+(3d0*n_children**2*A**2*B-1d0)*myroots**2+3*n_children**2*A*B**2*myroots+n_children**2*B**3
00110          !      stop
00111          IF (j == 0) THEN
00112             IF (MPI_ID == 0) WRITE(*,*) 'No Real Roots Found.  Just ', MyRoots
00113             STOP
00114          ELSEIF (j == 1) THEN
00115             IF (MPI_ID == 0) WRITE(*,*) 'Found 1 Real Root', ValidRoots(1)
00116          ELSE
00117             IF (MPI_ID == 0) WRITE(*,*) 'Found several real roots... using largest of ValidRoots(1:j)'
00118             ValidRoots(1)=minval(ValidRoots(1:j))
00119          END IF
00120          fill_fraction(n)=ValidRoots(1)      
00121          xi(n)=chi*xi(n+1)*fill_fraction(n)+(1d0-fill_fraction(n)) !over density of level n-1 due to clumplets
00122       END DO
00123       r_children(0)=r_cloud
00124       DO n=1, ClumpTreeDepth
00125          r_children(n)=r_children(n-1)*fill_fraction(n-1)**(1d0/3d0) !REAL(nDim))
00126       END DO
00127       DO n=0,ClumpTreeDepth
00128          rho_children(n)=JeansDensity(r_children(n), MinTemp)/xi(n)
00129       END DO
00130 
00131 !      write(*,*) rho_children(0)*(1d0-fill_fraction(0))+fill_fraction(0)*((1d0-fill_fraction(1))*chi+fill_fraction(1)*chi*chi)/xi(0)
00132       write(*,*) 'xi=', xi(0:ClumpTreeDepth)
00133       write(*,*) 'r_children=', r_children(0:ClumpTreeDepth)
00134       write(*,*) 'rho_children=', rho_children(0:ClumpTreeDepth)
00135       write(*,*) 'filling_fraction', fill_fraction(0:ClumpTreeDepth)
00136       CALL CreateAmbient(Ambient)
00137       Ambient%density=JeansDensity(r_ambient,MinTemp)
00138 
00139       CLOSE(PROBLEM_DATA_HANDLE)      
00140       
00141      if(lUseExisting) then    
00142         OPEN(UNIT=PROBLEM_DATA_HANDLE, FILE='clumplets.data', STATUS="OLD")
00143         write(*,*) sum(n_children**(/(n, n=0,ClumpTreeDepth)/))
00144         do i=1,sum(n_children**(/(n, n=0,ClumpTreeDepth)/))
00145            CALL CreateClump(Clump)
00146            read(PROBLEM_DATA_HANDLE,'(5e25.16)') Clump%density, Clump%position, Clump%radius
00147            CALL UpdateClump(Clump)
00148         end do 
00149      else
00150         OPEN(UNIT=PROBLEM_DATA_HANDLE, FILE='clumplets.data', STATUS="UNKNOWN")
00151         CALL CreateClumplets(half*Sum(GxBounds(:,:),2), 0) !puts it in middle
00152      end if
00153 
00154       CLOSE(PROBLEM_DATA_HANDLE)      
00155 
00156       DO i=1,nDim
00157          DO edge=1,2
00158             IF (Gmthbc(i,edge) == 1) THEN 
00159                CALL CreateWind(Wind)
00160                Wind%dir=i
00161                Wind%edge=edge
00162                Wind%type=OUTFLOW_ONLY
00163                Wind%density=Ambient%density
00164                Wind%temperature=Ambient%Pressure/Ambient%Density
00165                Wind%B=Ambient%B
00166             END IF
00167          END DO
00168       END DO
00169 
00170 
00171 !      CALL AddAllTotals(GASCOMP)
00172 !      CALL AddAllTotals(PARTICLECOMP)
00173 !      CALL AddAllTotals(BOTHCOMP)
00174 
00175 !      CALL CreateHistogram(Histogram)
00176 !      Histogram%Field%iD=1
00177 !      Histogram%Field%component=GASCOMP
00178 !      Histogram%minvalue=.1d0
00179 !      Histogram%maxvalue=1d8
00180       !      Histogram%nbins=nbins
00181 !      Histogram%scale=LOGSCALE
00182 
00183 !      CALL CreateProjection(projection)
00184 !      Projection%Field%iD=Mass_Field
00185 !      Projection%Field%component=BOTHCOMP
00186 !      Projection%dim=3
00187 !      ALLOCATE(Projection%Image)
00188 !      Projection%Image%Scaling=LOGSCALE
00189 
00190 !      CALL CreateProjection(projection)
00191 !      Projection%Field%iD=Mass_Field
00192 !      Projection%Field%component=BOTHCOMP
00193 !      Projection%dim=1
00194 !      ALLOCATE(Projection%Image)
00195 !      Projection%Image%Scaling=LOGSCALE
00196 
00197 !      CALL CreateProjection(projection)
00198 !      Projection%Field%iD=Mass_Field
00199 !      Projection%Field%component=BOTHCOMP
00200 !      Projection%pow=1d0
00201 !      Projection%dim=2
00202 !      ALLOCATE(Projection%Image)
00203 !      Projection%Image%Scaling=LOGSCALE
00204 
00205       CALL ClearAllRefinements()
00206       CALL AddRefinementThreshold(JeansLength_Field, LESSTHAN, (/(16d0*levels(i)%dx,i=0,MaxLevel)/))
00207 
00208    END SUBROUTINE ProblemModuleInit
00209 
00210    RECURSIVE SUBROUTINE CreateClumplets(pos,n)
00211       REAL(KIND=qPREC), DIMENSION(3) :: pos
00212       INTEGER :: n,i
00213       TYPE(ClumpDef), POINTER :: Clump
00214       REAL(KIND=qPREC), DIMENSION(:,:), ALLOCATABLE :: childpos
00215 
00216       CALL CreateClump(Clump)
00217       Clump%density=rho_children(n)
00218       Clump%position=pos
00219       Clump%radius=r_children(n)
00220       CALL UpdateClump(Clump)
00221 
00222       if (MPI_ID == 0) then
00223          write(PROBLEM_DATA_HANDLE,'(5e25.16)') Clump%density, Clump%position, Clump%radius
00224       end if
00225 
00226       IF (n < ClumpTreeDepth) THEN
00227          ALLOCATE(childpos(n_children,3))
00228          CALL GetChildPos(n_children,childpos,1.0*fill_fraction(n)**(1d0/3d0)) !REAL(nDim)))
00229          DO i=1,n_children
00230             CALL CreateClumpLets(pos+r_children(n)*childpos(i,:), n+1)
00231          END DO
00232          DEALLOCATE(childpos)
00233       END IF
00234    END SUBROUTINE CreateClumplets
00235 
00236 
00239   SUBROUTINE ProblemGridInit(Info)
00240     TYPE(InfoDef) :: Info
00241   END SUBROUTINE ProblemGridInit
00242 
00245   SUBROUTINE ProblemBeforeStep(Info)
00246     TYPE(InfoDef) :: Info
00247  END SUBROUTINE ProblemBeforeStep
00248 
00251   SUBROUTINE ProblemAfterStep(Info)
00252     TYPE(InfoDef) :: Info
00253 
00254   END SUBROUTINE ProblemAfterStep
00255 
00256 
00259   SUBROUTINE ProblemSetErrFlag(Info)
00260     USE CoolingSrc
00261     TYPE(InfoDef) :: Info
00262     
00263  END SUBROUTINE ProblemSetErrFlag
00264 
00265 
00266   SUBROUTINE ProblemBeforeGlobalStep(n)
00267      INTEGER :: n
00268   END SUBROUTINE ProblemBeforeGlobalStep
00269 
00270 
00271   FUNCTION GetRoots(a,b,c,d,nroot)
00272      COMPLEX(8) :: GetRoots(3)
00273      REAL(8) :: a,b,c,d,DD,p,q,u,v,temp1,temp2,phi
00274      REAL(8) :: y1,y2,y3,y2r,y2i
00275      INTEGER :: nroot
00276      ! Variables used:
00277      !   a, b, c, d  ... coefficients (input)
00278      !   y1, y2, y3  ... three transformed solutions
00279      !   y2r, y2i    ... real and imaginary parts of a pair of complex roots
00280      !   nroot       ... number of roots
00281      !
00282      ! Formula used are given in Tuma, "Engineering Mathematics Handbook", p7
00283      !   (McGraw Hill, 1978).
00284      !   Step 0: If a is 0. use the quadrati! formula to avoid dividing by 0.
00285      !   Step 1: Calculate p and q
00286      !           p = ( 3*c/a - (b/a)**2 ) / 3
00287      !           q = ( 2*(b/a)**3 - 9*b*c/a/a + 27*d/a ) / 27
00288      !   Step 2: Calculate discriminant D
00289      !           D = (p/3)**3 + (q/2)**2
00290      !   Step 3: Depending on the sign of D, we follow different strategy.
00291      !           If D<0, three distinct real roots.
00292      !           If D=0, three real roots of which at least two are equal.
00293      !           If D>0, one real and two complex roots.
00294      !   Step 3a: For D>0 and D=0,
00295      !           Calculate u and v
00296      !           u = cubic_root(-q/2 + sqrt(D))
00297      !           v = cubic_root(-q/2 - sqrt(D))
00298      !           Find the three transformed roots
00299      !           y1 = u + v
00300      !           y2 = -(u+v)/2 + i (u-v)*sqrt(3)/2
00301      !           y3 = -(u+v)/2 - i (u-v)*sqrt(3)/2
00302      !   Step 3b Alternately, for D<0, a trigonometri! formulation is more convenient
00303      !           y1 =  2 * sqrt(|p|/3) * cos(phi/3)
00304      !           y2 = -2 * sqrt(|p|/3) * cos((phi+pi)/3)
00305      !           y3 = -2 * sqrt(|p|/3) * cos((phi-pi)/3)
00306      !           where phi = acos(-q/2/sqrt(|p|**3/27))
00307      !                 pi  = 3.141592654...
00308      !   Step 4  Finally, find the three roots
00309      !           x = y - b/a/3
00310      !
00311      ! ----------------------------------------------------------------------
00312      ! Instructor: Nam Sun Wang
00313      ! ----------------------------------------------------------------------
00314 
00315      ! Declare variables
00316 
00317      ! Step 0: If a is 0 use the quadrati! formula. -------------------------
00318      IF(a .eq. 0d0)THEN
00319         if(b .eq. 0d0)then
00320            if(c .eq. 0d0)then
00321               !         We have a non-equation; therefore, we have no valid solution
00322               nroot = 0
00323            else
00324               !         We have a linear equation with 1 root.
00325               nroot = 1
00326               Getroots(1) = cmplx(-d/c, 0d0)
00327            endif
00328         else
00329            !     We have a true quadrati! equation.  Apply the quadrati! formula to find two roots.
00330            nroot = 2
00331            DD = c*c-4d0*b*d
00332            if(DD .ge. 0d0)then
00333               Getroots(1) = cmplx((-c+sqrt(DD))/2d0/b, 0d0)
00334               Getroots(2) = cmplx((-c-sqrt(DD))/2d0/b, 0d0)
00335            else
00336               Getroots(1) = cmplx(-c/2d0/b, +sqrt(-DD)/2d0/b)
00337               Getroots(2) = cmplx(-c/2d0/b, -sqrt(-DD)/2d0/b)
00338            endif
00339         endif
00340 
00341      ELSE
00342 
00343         ! Cubi! equation with 3 roots
00344         nroot = 3
00345 
00346         ! Step 1: Calculate p and q --------------------------------------------
00347         p  = c/a - b*b/a/a/3d0
00348         q  = (2d0*b*b*b/a/a/a - 9d0*b*c/a/a + 27d0*d/a) / 27d0
00349 
00350         ! Step 2: Calculate DD (discriminant) ----------------------------------
00351         DD = p*p*p/27d0 + q*q/4d0
00352 
00353         ! Step 3: Branch to different algorithms based on DD ------------------
00354 
00355         if(DD .lt. 0d0)then
00356            !       Step 3b:
00357            !       3 real unequal roots -- use the trigonometri! formulation
00358            phi = acos(-q/2d0/sqrt(abs(p*p*p)/27d0))
00359            temp1=2d0*sqrt(abs(p)/3d0)
00360            y1 =  temp1*cos(phi/3d0)
00361            y2 = -temp1*cos((phi+pi)/3d0)
00362            y3 = -temp1*cos((phi-pi)/3d0)
00363            temp1 = b/a/3d0
00364            y2 = y2-temp1
00365            y3 = y3-temp1
00366 
00367         else
00368            !       Step 3a:
00369            !       1 real root & 2 conjugate complex roots OR 3 real roots (some are equal)
00370            temp1 = -q/2d0 + sqrt(DD)
00371            temp2 = -q/2d0 - sqrt(DD)
00372            u = abs(temp1)**(1d0/3d0)
00373            v = abs(temp2)**(1d0/3d0)
00374            if(temp1 .lt. 0d0) u=-u
00375            if(temp2 .lt. 0d0) v=-v
00376            y1  = u + v
00377            y2r = -(u+v)/2d0
00378            y2i =  (u-v)*sqrt(3d0)/2d0
00379            temp1 = b/a/3d0
00380            y2r=y2r-temp1
00381 
00382 
00383         endif
00384 
00385         ! Step 4: Final transformation -----------------------------------------
00386 
00387         y1 = y1-temp1
00388 
00389         ! Assign answers -------------------------------------------------------
00390         if(DD .lt. 0d0)then
00391            Getroots(1) = cmplx( y1,  0d0)
00392            Getroots(2) = cmplx( y2,  0d0)
00393            Getroots(3) = cmplx( y3,  0d0)
00394         elseif(DD .eq. 0d0)then
00395            Getroots(1) = cmplx( y1,  0d0)
00396            Getroots(2) = cmplx(y2r,  0d0)
00397            Getroots(3) = cmplx(y2r,  0d0)
00398         else
00399            Getroots(1) = cmplx( y1,  0d0)
00400            Getroots(2) = cmplx(y2r, y2i)
00401            Getroots(3) = cmplx(y2r,-y2i)
00402         endif
00403 
00404      ENDIF
00405 
00406   end FUNCTION GetRoots
00407 
00408 
00409 
00410   FUNCTION JeansDensity(r, temp)
00411      REAL(KIND=qPREC) :: r, temp, JeansDensity
00412      JeansDensity=gamma*temp/r**2*pi/ScaleGrav
00413   END FUNCTION JeansDensity
00414 
00415 
00416   SUBROUTINE GetChildPos(n_children,childpos,r)
00417      INTEGER :: n_children
00418      REAL(KIND=qPREC), DIMENSION(:,:) :: childpos
00419      REAL(KIND=qPREC) :: r
00420      INTEGER, PARAMETER :: nIterations=10000
00421      LOGICAL :: lSuccess
00422      INTEGER :: j,k,l
00423      REAL :: a
00424      !create childpositions where each child is at least 2r apart and within (1-r) of the center...
00425      childpos=0
00426      DO j=1,nIterations
00427         lSuccess=.true.
00428         DO k=1,n_children
00429            CALL Random(a)
00430            childpos(k,1)=a
00431            CALL Random(a)
00432            childpos(k,2)=a
00433            IF (nDim == 3) THEN
00434               CALL Random(a)
00435               childpos(k,3)=a
00436            END IF
00437 !           write(*,*) 'try', Childpos(k,:)
00438            Childpos(k,1:nDim)=2d0*(childpos(k,1:nDim)-.5d0)*(1d0-r) !rescale to be within a random cube with half width (1-r)                      
00439 !           write(*,*) 'trying', Childpos(k,:)
00440            IF (sqrt(sum(childpos(k,1:nDim)**2)) > (1d0-r)) THEN
00441 !              write(*,*) 'too close to boudnary'
00442               lSuccess=.false.
00443            ELSE
00444               DO l=1,k-1
00445                  IF (sqrt(sum((childpos(k,1:nDim)-childpos(l,1:nDim))**2)) < 2d0*r) THEN
00446                     lSuccess=.false.
00447 !                    write(*,*) 'too close to sibling'
00448                     EXIT
00449                  END IF
00450               END DO
00451            END IF
00452            IF (.Not. lSuccess) EXIT
00453         END DO
00454         IF (lSuccess) EXIT
00455      END DO
00456      IF (.NOT. lSuccess) THEN
00457         IF (MPI_ID == 0) write(*,*) 'failed to place', n_children, 'clumplets of radius ',r
00458         
00459         STOP
00460      END IF
00461 
00462   END SUBROUTINE GetChildPos
00463 
00464 
00465 
00466 
00467 END MODULE Problem
00468 
 All Classes Files Functions Variables