Scrambler
1
|
00001 !######################################################################### 00002 ! 00003 ! Copyright (C) 2003-2012 Department of Physics and Astronomy, 00004 ! University of Rochester, 00005 ! Rochester, NY 00006 ! 00007 ! shapes.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 00029 00032 MODULE Shapes 00033 USE GlobalDeclarations 00034 USE PhysicsDeclarations 00035 USE ObjectDeclarations 00036 IMPLICIT NONE 00037 INTEGER, PARAMETER :: SPHERE=1, ELLIPSOID=2, CYLINDER=3, ELLIPTICAL_PRISM=4, RECTANGULAR_PRISM=5 00038 00039 TYPE ShapeDef 00040 INTEGER :: Type = SPHERE 00041 REAL(KIND=qPREC) :: position(3)=(/0d0,0d0,0d0/) 00042 ! REAL(KIND=qPREC) :: psi=0 !Initial rotation about z axis 00043 ! REAL(KIND=qPREC) :: theta=0 !Angle rotated around y (from z towards x) 00044 ! REAL(KIND=qPREC) :: phi=0 !Angle rotated around z (from x towards y) 00045 INTEGER :: xy(3) = (/1,2,3/) !integer indices corresponding to frame reference x and y directions 00046 REAL(KIND=qPREC) :: size_param(3) = (/1d0,1d0,1d0/) 00047 REAL(KIND=qPREC) :: RotationMatrix(3,3)=RESHAPE((/1d0,0d0,0d0,0d0,1d0,0d0,0d0,0d0,1d0/),(/3,3/)) !Rotates lab frame into object frame 00048 REAL(KIND=qPREC) :: xBounds(3,2)=RESHAPE((/-1d0,-1d0,-1d0,1d0,1d0,1d0/),(/3,2/)) 00049 REAL(KIND=qPREC) :: velocity(3)=0 00050 REAL(KIND=qPREC) :: t0 = 0 00051 LOGICAL :: offcenter = .false. 00052 REAL(KIND=qPREC) :: bounds_param(3,2) = RESHAPE((/-1d0,-1d0,-1d0,1d0,1d0,1d0/),(/3,2/)) 00053 00054 END TYPE ShapeDef 00055 00056 CONTAINS 00057 00058 SUBROUTINE CreateShape(Shape, type, dimensions) 00059 TYPE(ShapeDef), POINTER :: Shape 00060 INTEGER, OPTIONAL :: type 00061 REAL(KIND=qPREC), OPTIONAL :: dimensions(3) 00062 ALLOCATE(Shape) 00063 IF(Present(type) .AND. Present(dimensions)) THEN 00064 CALL SetShapeType(shape, type, dimensions) 00065 ENDIF 00066 CALL SetShapeBounds(Shape) 00067 END SUBROUTINE CreateShape 00068 00069 SUBROUTINE UpdateShape(Shape) 00070 TYPE(ShapeDef), POINTER :: Shape 00071 CALL SetShapeBounds(Shape) 00072 END SUBROUTINE UpdateShape 00073 00074 SUBROUTINE SetShapeOrientation(Shape, psi, theta, phi) 00075 TYPE(ShapeDef), POINTER :: Shape 00076 REAL(KIND=qPREC) :: c1,c2,c3,s1,s2,s3 00077 REAL(KIND=qPREC) :: theta,phi,psi 00078 REAL(KIND=qPREC) :: theta_,phi_,psi_ 00079 00080 00081 IF (nDim == 2) THEN 00082 !theta must be 0 Pi/2 or Pi 00083 !psi should be 0 Pi/2 Pi, 3Pi/2, 2Pi, etc... 00084 !phi should be 0 Pi/2 Pi, etc... 00085 ! psi_=Pi/2d0*nint(2d0*psi/Pi) 00086 ! theta_=Pi/2d0*nint(2d0*theta/Pi) 00087 ! phi_=Pi/2d0*nint(2d0*phi/Pi) 00088 ! write(*,*) 'adjusted psi, theta, and phi to', psi, theta, phi 00089 00090 psi_=psi 00091 theta_=theta 00092 phi_=phi 00093 00094 ELSE 00095 psi_=psi 00096 theta_=theta 00097 phi_=phi 00098 END IF 00099 00100 c1=cos(psi_) 00101 c2=cos(theta_) 00102 c3=cos(phi_) 00103 s1=sin(psi_) 00104 s2=sin(theta_) 00105 s3=sin(phi_) 00106 00107 00108 IF(abs(real(nint(c1),8)-c1)<1d-8) c1=real(nint(c1),8) 00109 IF(abs(real(nint(c2),8)-c2)<1d-8) c2=real(nint(c2),8) 00110 IF(abs(real(nint(c3),8)-c3)<1d-8) c3=real(nint(c3),8) 00111 IF(abs(real(nint(s1),8)-s1)<1d-8) s1=real(nint(s1),8) 00112 IF(abs(real(nint(s2),8)-s2)<1d-8) s2=real(nint(s2),8) 00113 IF(abs(real(nint(s3),8)-s3)<1d-8) s3=real(nint(s3),8) 00114 ! Shape%RotationMatrix=RESHAPE((/ & 00115 ! c1*c2*c3-s1*s3, -c2*s3*c1-c3*s1, c1*s2,& 00116 ! c1*s3+c3*c2*s1, c1*c3-c2*s1*s3, s2*s1, & 00117 ! -c3*s2, s3*s2, c2/),(/3,3/),(/0d0/),(/2,1/)) 00118 00119 00120 Shape%RotationMatrix=RESHAPE((/ & 00121 c1*c2*c3-s1*s3, -c3*c2*s1-c1*s3, +c3*s2,& 00122 c1*c2*s3+c3*s1, c1*c3-c2*s1*s3, s2*s3, & 00123 -c1*s2, s1*s2, c2/),(/3,3/),(/0d0/),(/2,1/)) 00124 00125 00126 ! IF (nDim == 2) THEN 00127 shape%xy(1)=maxval(maxloc(abs(RotateVectorToShape(shape,(/1d0,0d0,0d0/))))) 00128 shape%xy(2)=maxval(maxloc(abs(RotateVectorToShape(shape,(/0d0,1d0,0d0/))))) 00129 shape%xy(3)=maxval(maxloc(abs(RotateVectorToShape(shape,(/0d0,0d0,1d0/))))) 00130 ! IF (MPI_ID == 0) write(*,*) 'shape%xy=', shape%xy 00131 ! END IF 00132 00133 CALL SetShapeBounds(Shape) 00134 END SUBROUTINE SetShapeOrientation 00135 00136 00137 SUBROUTINE SetShapeType(Shape, type, dimensions) 00138 TYPE(ShapeDef), POINTER :: Shape 00139 INTEGER :: type 00140 REAL(KIND=qPREC) :: dimensions(3) 00141 Shape%type=type 00142 SELECT CASE (type) 00143 CASE (SPHERE) 00144 Shape%size_param=dimensions(1) 00145 CASE (CYLINDER) 00146 Shape%size_param=(/dimensions(1),dimensions(1),dimensions(3)/) 00147 CASE (ELLIPSOID, ELLIPTICAL_PRISM, RECTANGULAR_PRISM) 00148 SHape%size_param=dimensions 00149 CASE DEFAULT 00150 PRINT*, 'unrecognized shape type in shapes.f90' 00151 STOP 00152 END SELECT 00153 END SUBROUTINE SetShapeType 00154 00155 FUNCTION TransformCoordsFromShape(Shape, pos, tnow) 00156 TYPE(ShapeDef), POINTER :: Shape 00157 REAL(KIND=qPREC), DIMENSION(3) :: TransformCoordsFromShape, pos 00158 REAL(KIND=qPREC), OPTIONAL :: tnow 00159 INTEGER :: l 00160 IF (PRESENT(tnow)) THEN 00161 DO l=1,3 00162 TransformCoordsFromShape(l)=Shape%position(l)+Shape%velocity(l)*(tnow-Shape%t0)+SUM(Shape%RotationMatrix(l,:)*(pos(:))) !Inverse rotation 00163 END DO 00164 00165 ELSE 00166 DO l=1,3 00167 TransformCoordsFromShape(l)=Shape%position(l)+SUM(Shape%RotationMatrix(l,:)*(pos(:))) !Inverse rotation 00168 END DO 00169 END IF 00170 END FUNCTION TransformCoordsFromShape 00171 00172 FUNCTION TransformCoordsToShape(Shape, pos, tnow) 00173 TYPE(ShapeDef), POINTER :: Shape 00174 REAL(KIND=qPREC), DIMENSION(3) :: TransformCoordsToShape, pos 00175 REAL(KIND=qPREC), OPTIONAL :: tnow 00176 INTEGER :: l 00177 IF (Present(tnow)) THEN 00178 DO l=1,3 00179 TransformCoordsToShape(l)=SUM(Shape%RotationMatrix(:,l)*(pos(:)-(Shape%position(:)+Shape%velocity*(tnow-Shape%t0)))) !Inverse rotation 00180 END DO 00181 ELSE 00182 DO l=1,3 00183 TransformCoordsToShape(l)=SUM(Shape%RotationMatrix(:,l)*(pos(:)-Shape%position(:))) !Inverse rotation 00184 END DO 00185 00186 END IF 00187 END FUNCTION TransformCoordsToShape 00188 00189 FUNCTION RotateVectorFromShape(Shape, pos) 00190 TYPE(ShapeDef), POINTER :: Shape 00191 REAL(KIND=qPREC), DIMENSION(3) :: RotateVectorFromShape, pos 00192 INTEGER :: l 00193 DO l=1,3 00194 RotateVectorFromShape(l)=SUM(Shape%RotationMatrix(l,:)*pos(:)) !Inverse rotation 00195 END DO 00196 END FUNCTION RotateVectorFromShape 00197 00198 FUNCTION RotateVectorToShape(Shape, pos) 00199 TYPE(ShapeDef), POINTER :: Shape 00200 REAL(KIND=qPREC), DIMENSION(3) :: RotateVectorToShape, pos 00201 INTEGER :: l 00202 DO l=1,3 00203 RotateVectorToShape(l)=SUM(Shape%RotationMatrix(:,l)*pos(:)) !Inverse rotation 00204 END DO 00205 END FUNCTION RotateVectorToShape 00206 00207 00208 SUBROUTINE SetShapeBounds(Shape) 00209 TYPE(ShapeDef), POINTER :: Shape 00210 REAL(KIND=qPREC), DIMENSION(3) :: pos !Position of each corner in object frame 00211 REAL(KIND=qPREC), DIMENSION(3) :: rpos !Rotated position in simulation frame 00212 INTEGER :: i,j,k 00213 Shape%xBounds=SPREAD(Shape%position,2,2) 00214 SELECT CASE (Shape%Type) 00215 CASE(SPHERE) !Orientation does not matter 00216 Shape%xBounds(:,1)=Shape%position-Shape%size_param(1) 00217 Shape%xBounds(:,2)=Shape%position+Shape%size_param(1) 00218 CASE(ELLIPSOID, CYLINDER, ELLIPTICAL_PRISM, RECTANGULAR_PRISM) !Bound shape by rectangle and find max extents of each corner. 00219 IF(.NOT. Shape%offcenter) THEN 00220 DO i=1,2 00221 pos(1)=(-1d0)**i*Shape%size_param(1) 00222 DO j=1,2 00223 pos(2)=(-1d0)**j*Shape%size_param(2) 00224 DO k=1,2 00225 pos(3)=(-1d0)**k*Shape%size_param(3) 00226 rpos=TransformCoordsFromShape(Shape, pos) 00227 Shape%xBounds(:,1)=min(Shape%xBounds(:,1),rpos) 00228 Shape%xBounds(:,2)=max(Shape%xBounds(:,2),rpos) 00229 END DO 00230 END DO 00231 END DO 00232 ELSE 00233 DO i=1,2 00234 pos(1) = Shape%bounds_param(1,i) 00235 DO j=1,2 00236 pos(2) = Shape%bounds_param(2,j) 00237 DO k=1,2 00238 pos(3) = Shape%bounds_param(3,k) 00239 rpos = TransformCoordsFromShape(Shape,pos) 00240 Shape%xBounds(:,1)=min(Shape%xBounds(:,1),rpos) 00241 Shape%xBounds(:,2)=max(Shape%xBounds(:,2),rpos) 00242 END DO 00243 END DO 00244 END DO 00245 END IF 00246 END SELECT 00247 END SUBROUTINE SetShapeBounds 00248 00249 FUNCTION GetShapeBounds(Shape, tnow) 00250 TYPE(ShapeDef) :: Shape 00251 REAL(KIND=qPREC), DIMENSION(3,2) :: GetShapeBounds 00252 REAL(KIND=qPREC), OPTIONAL :: tnow 00253 IF (Present(tnow)) THEN 00254 GetShapeBounds=Shape%xBounds+SPREAD(Shape%velocity*(tnow-Shape%t0),2,2) 00255 ELSE 00256 GetShapeBounds=Shape%xBounds 00257 END IF 00258 END FUNCTION GetShapeBounds 00259 00260 FUNCTION IsInShape(Shape, pos, coords, tnow) 00261 TYPE(ShapeDef), POINTER :: Shape 00262 REAL(KIND=qPREC), DIMENSION(3) :: pos 00263 REAL(KIND=qPREC), DIMENSION(3) :: coords 00264 REAL(KIND=qPREC), DIMENSION(3) :: rpos 00265 REAL(KIND=qPREC), OPTIONAL :: tnow 00266 LOGICAL :: IsInShape 00267 IsInShape=.false. 00268 IF (PRESENT(tnow)) THEN 00269 coords=TransformCoordsToShape(Shape, pos, tnow) 00270 ELSE 00271 coords=TransformCoordsToShape(Shape, pos) 00272 END IF 00273 SELECT CASE (SHAPE%TYPE) 00274 CASE(SPHERE) 00275 IsInShape=(sum(coords(:)**2) <= Shape%size_param(1)**2) 00276 CASE (CYLINDER) 00277 IsInShape=(sum(coords(1:2)**2) <= Shape%size_param(1)**2) .AND. (abs(coords(3)) <= Shape%size_param(3)) 00278 CASE (RECTANGULAR_PRISM) 00279 IsInShape=ALL(abs(coords(:)) <= Shape%size_param(:)) 00280 CASE (ELLIPSOID) 00281 IsInShape=(SUM(coords(:)**2/Shape%size_param(:)**2) <= 1d0) 00282 CASE (ELLIPTICAL_PRISM) 00283 IsInShape=(SUM(coords(1:2)**2/Shape%size_param(1:2)**2) <= 1d0) .AND. (abs(coords(3)) <= Shape%size_param(3)) 00284 END SELECT 00285 00286 END FUNCTION IsInShape 00287 00288 SUBROUTINE PrintShape(shape, tnow) 00289 TYPE(ShapeDef), POINTER :: Shape 00290 REAL(KIND=qPREC), OPTIONAL :: tnow 00291 REAL(KIND=qPREC) :: pos(3) 00292 IF (PRESENT(tnow)) THEN 00293 pos=Shape%position+Shape%velocity*(tnow-shape%t0) 00294 ELSE 00295 pos=Shape%position 00296 END IF 00297 write(*,*) '================= Shape ====================' 00298 SELECT CASE (shape%Type) 00299 CASE(SPHERE) 00300 write(*,'(A,E24.15,A,3E24.15)') 'Sphere of radius ', shape%size_param(1), 'centered at ', pos(1:nDim) 00301 CASE(CYLINDER) 00302 write(*,'(A,E24.15,A,E24.15,A,3E24.15)') 'Cylinder of radius', shape%size_param(1), ' and height ', 2d0*shape%size_param(3), ' centered at ', pos(1:nDim), ' with axial vector ', RotateVectorFromShape(shape,(/0d0,0d0,01d0/)) 00303 CASE(RECTANGULAR_PRISM) 00304 write(*,'(A,E24.15,A,E24.15,A,E24.15,A,3E24.15)') 'Rectangular prism of width', 2d0*shape%size_param(1), ' height ', 2d0*shape%size_param(3), ' and depth ', 2d0*shape%size_param(2), ' centered at ', pos(1:nDim) 00305 write(*,'(A)') 'Vectors normal to the x, y, and z faces are: ' 00306 write(*,'(3E24.15)') RotateVectorFromShape(shape,(/1d0,0d0,0d0/)), RotateVectorFromShape(shape,(/0d0,1d0,0d0/)), RotateVectorFromShape(shape,(/0d0,0d0,1d0/)) 00307 CASE(ELLIPSOID) 00308 write(*,'(A,E24.15,A,E24.15,A,E24.15,A,3E24.15)') 'Ellipsoid of width', 2d0*shape%size_param(1), ' height ', 2d0*shape%size_param(3), ' and depth ', 2d0*shape%size_param(2), ' centered at ', pos(1:nDim) 00309 write(*,'(A)') 'Vectors corresponding to the x, y, and z axis are: ' 00310 write(*,'(3E24.15)') RotateVectorFromShape(shape,(/1d0,0d0,0d0/)), RotateVectorFromShape(shape,(/0d0,1d0,0d0/)), RotateVectorFromShape(shape,(/0d0,0d0,1d0/)) 00311 CASE(ELLIPTICAL_PRISM) 00312 write(*,'(A,2E24.15,A,E24.15,A,3E24.15)') 'Elliptical prism with radii',shape%size_param(1:2), ' and height ', 2d0*shape%size_param(3), ' centered at ', pos(1:nDim) 00313 write(*,'(A)') 'Vectors corresponding to the x, y, and z axis are: ' 00314 write(*,'(3E24.15)') RotateVectorFromShape(shape,(/1d0,0d0,0d0/)), RotateVectorFromShape(shape,(/0d0,1d0,0d0/)), RotateVectorFromShape(shape,(/0d0,0d0,1d0/)) 00315 00316 END SELECT 00317 write(*,*) '============================================' 00318 END SUBROUTINE PrintShape 00319 00320 SUBROUTINE GetRandomPointsInShape(shape, min_separation, boundary_separation, positions, tnow) 00321 TYPE(ShapeDef), POINTER :: shape 00322 INTEGER :: npoints, recheck_count, j,k,l 00323 REAL(KIND=qPREC) :: min_separation,sep,minsep,boundary_separation 00324 REAL :: a 00325 REAL(KIND=qPREC), DIMENSION(:,:) :: positions 00326 REAL(KIND=qPREC), OPTIONAL :: tnow !current time (needed if shape is moving) 00327 INTEGER, PARAMETER :: MaxIter=1000, MaxRetries=10 00328 LOGICAL :: lSuccess, lReCheck 00329 INTEGER :: ii,jj,kk,ioffset(3,2), iErr 00330 REAL(KIND=qPREC), DIMENSION(3) :: poffset, minpos,temp,temp2, periodic_pos 00331 ioffset=0 00332 WHERE (lHydroPeriodic(1:nDim)) ioffset(1:nDim,2)=1 00333 ioffset(1:nDim,1)=-ioffset(1:nDim,2) 00334 00335 ! write(*,*) min_separation, boundary_separation, shape%size_param 00336 npoints=size(positions,1) 00337 positions=0d0 !spread(shape%position, 1, npoints) 00338 DO k=1,npoints 00339 DO j=1,MaxIter 00340 lSuccess=.true. 00341 CALL Random_number(a) 00342 positions(k,shape%xy(1))=(a*2d0-1d0)*(shape%size_param(shape%xy(1))-boundary_separation) 00343 CALL Random_number(a) 00344 positions(k,shape%xy(2))=(a*2d0-1d0)*(shape%size_param(shape%xy(2))-boundary_separation) 00345 IF (nDim == 3) THEN 00346 CALL Random_number(a) 00347 positions(k,shape%xy(3))=(a*2d0-1d0)*(shape%size_param(shape%xy(3))-boundary_separation) 00348 END IF 00349 00350 IF (MPI_ID == 0) writE(*,*) 'placing particle ', k, 'of ', npoints ![Bpositions(k,:) 00351 recheck_count=0 !how many times have we tried shifting the point away from the boundary and other points 00352 lReCheck=.true. !Do we need to recheck the distance from other points after moving this one? 00353 DO WHILE (lReCheck .AND. recheck_count < MaxReTries) 00354 lSuccess=.true. !Assume success in finding a point with no conflicts 00355 lReCheck=.false. !Assume that we will not shift the point 00356 DO l=1,k-1 ! Loop over all current existing points 00357 minsep=1d30 00358 minpos=positions(l,:) 00359 DO ii=ioffset(1,1),ioffset(1,2) 00360 DO jj=ioffset(2,1),ioffset(2,2) 00361 DO kk=ioffset(3,1),ioffset(3,2) 00362 pOffSet=(/ii,jj,kk/)*(GxBounds(:,2)-GxBounds(:,1)) 00363 periodic_pos=positions(l,:)+RotateVectorToShape(shape,pOffSet) 00364 sep=sqrt(sum((positions(k,:)-periodic_pos)**2)) !distance from point to existing point 00365 IF (sep < minsep) THEN 00366 ! IF (MPI_ID == 0) write(*,'(A1,2I,7E24.15)') 'a',k,l,sep,positions(k,:),periodic_pos 00367 minpos=periodic_pos 00368 minsep=sep 00369 END IF 00370 END DO 00371 END DO 00372 END DO 00373 temp=positions(k,:) 00374 ! IF (MPI_ID == 0) write(*,*) 'temp=', temp 00375 IF (minsep < min_separation) THEN 00376 !try shifting it away from that particle so it is 1.1 min_separations away 00377 ! make sure it is still inside the shape bounds 00378 ! IF (MPI_ID == 0) write(*,*) 'shifting from particle', l, positions(l,:), minpos, minsep 00379 positions(k,shape%xy(1:nDim))=positions(k,shape%xy(1:nDim))+min_separation/minsep*1.1d0*(positions(k,shape%xy(1:nDim))-minpos(shape%xy(1:nDim))) 00380 END IF 00381 temp2=positions(k,:) 00382 ! IF (MPI_ID == 0) write(*,*) 'temp2=', temp2 00383 IF (ANY(temp /= positions(k,:)) .OR. l == 1) THEN 00384 SELECT CASE(Shape%Type) 00385 CASE(SPHERE) 00386 positions(k,shape%xy(1:nDim))=positions(k,shape%xy(1:nDim))*min(1d0,(shape%size_param(1)-boundary_separation)/sqrt(sum(positions(k,shape%xy(1:nDim))**2))) 00387 CASE (CYLINDER) 00388 IF (nDim == 2) THEN 00389 IF (shape%xy(3) == 3) THEN !cylinder projects to a circle in 2D 00390 positions(k,shape%xy(1:2))=positions(k,shape%xy(1:2))*min(1d0,(shape%size_param(1)-boundary_separation)/sqrt(sum(positions(k,shape%xy(1:2))**2))) 00391 ELSE 00392 positions(k,shape%xy(1:2))=positions(k,shape%xy(1:2))*min(1d0, (shape%size_param(shape%xy(1:2))-boundary_separation)/abs(positions(k,shape%xy(1:2)))) 00393 END IF 00394 ELSE 00395 positions(k,shape%xy(1:2))=positions(k,shape%xy(1:2))*min(1d0, (shape%size_param(1)-boundary_separation)/sqrt(sum(positions(k,shape%xy(1:2))**2))) 00396 positions(k,shape%xy(3))=positions(k,shape%xy(3))*min(1d0, (shape%size_param(3)-boundary_separation)/abs(positions(k,shape%xy(3)))) 00397 END IF 00398 CASE (RECTANGULAR_PRISM) 00399 positions(k,shape%xy(1:nDim))=positions(k,shape%xy(1:nDim))*min(1d0,(shape%size_param(shape%xy(1:nDim))-boundary_separation)/abs(positions(k,shape%xy(1:nDim)))) 00400 CASE (ELLIPSOID, ELLIPTICAL_PRISM) 00401 IF (Shape%Type == ELLIPSOID .AND. ALL(Shape%size_param(shape%xy(1:nDim))==maxval(shape%size_param(shape%xy(1:nDim))))) THEN 00402 positions(k,shape%xy(1:nDim))=positions(k,shape%xy(1:nDim))*min(1d0,(shape%size_param(1)-boundary_separation)/sqrt(sum(positions(k,shape%xy(1:nDim)**2)))) 00403 ELSEIF (Shape%Type == ELLIPTICAL_PRISM .AND. (nDim == 2 .OR. ALL(Shape%size_param(1:2)==maxval(shape%size_param(1:2))))) THEN 00404 IF (nDim == 2) THEN 00405 IF (shape%xy(3) == 3) THEN 00406 positions(k,1:2)=positions(k,1:2)*min(1d0,(shape%size_param(1)-boundary_separation)/sqrt(sum(positions(k,1:2)**2))) 00407 ELSE 00408 positions(k,shape%xy(1:2))=positions(k,shape%xy(1:2))*min(1d0,(shape%size_param(shape%xy(1:2))-boundary_separation)/abs(positions(k,shape%xy(1:2)))) 00409 END IF 00410 ELSE 00411 positions(k,1:2)=positions(k,1:2)*min(1d0,(shape%size_param(1)-boundary_separation)/sqrt(sum(positions(k,1:2)**2))) 00412 positions(k,3)=positions(k,3)*min(1d0,(shape%size_param(3)-boundary_separation)/abs(positions(k,3))) 00413 END IF 00414 ELSE 00415 IF (MPI_ID == 0) THEN 00416 write(*,*) 'getrandompointsinshape not implemented for ellipsoids or elliptical prisms yet that are not effectively spheres/circles/ellipses/cylinders. ' 00417 write(*,*) 'Prescription is available online at http://www.geometrictools.com/Documentation/DistancePointEllipseEllipsoid.pdf ' 00418 write(*,*) 'STOPPING' 00419 END IF 00420 CALL MPI_BARRIER(MPI_COMM_WORLD, iErr) 00421 STOP 00422 END IF 00423 END SELECT 00424 END IF 00425 IF (ANY(temp2 /= positions(k,:))) THEN 00426 ! IF (MPI_ID == 0) write(*,*) 'moved away from boundary' 00427 END IF 00428 IF (ANY(temp /= positions(k,:))) THEN 00429 ! IF (MPI_ID == 0) write(*,*) 'rechecking', temp, temp2, positions(k,:) 00430 lRecheck=.true. 00431 recheck_count=recheck_count+1 00432 lSuccess=.false. 00433 EXIT 00434 END IF 00435 END DO 00436 !Do we need to recheck? 00437 END DO 00438 ! Have we succeeded in this particle 00439 IF (lSuccess) EXIT 00440 END DO 00441 IF (.NOT. lSuccess) THEN 00442 IF (MPI_ID == 0) THEN 00443 write(*,'(A,I6,A,I6,A,E24.15, A)') 'Only was able to place ', k-1 , ' of ', npoints, ' points with a min separation of ', min_separation, ' inside of : ' 00444 CALL PrintShape(shape) 00445 write(*,'(A,E24.15,A)') 'min_separation is ', min_separation/(shapevolume(shape)/npoints)**(1d0/REAL(nDim)), ' of the grid packing value. Try reducing the number of points or the min_separation' 00446 END IF 00447 CALL MPI_BARRIER(MPI_COMM_WORLD, iErr) 00448 STOP 00449 END IF 00450 END DO 00451 DO k=1,npoints 00452 IF (PRESENT(tnow)) THEN 00453 positions(k,:)=TransformCoordsFromShape(shape, positions(k,:), tnow) 00454 ELSE 00455 positions(k,:)=TransformCoordsFromShape(shape, positions(k,:)) 00456 END IF 00457 END DO 00458 ! write(*,*) positions(:,:) 00459 ! STOP 00460 00461 END SUBROUTINE GetRandomPointsInShape 00462 00463 00464 FUNCTION ShapeVolume(shape) 00465 REAL(KIND=qPREC) :: ShapeVolume 00466 TYPE(ShapeDef), POINTER :: shape 00467 SELECT CASE(Shape%TYPE) 00468 CASE(SPHERE) 00469 IF (nDim == 2) THEN 00470 ShapeVolume=pi*shape%size_param(1)**2 00471 ELSE 00472 ShapeVolume=4d0/3d0*pi*shape%size_param(1)**3 00473 END IF 00474 CASE (CYLINDER) 00475 IF (nDim == 2) THEN 00476 IF (shape%xy(3) == 3) THEN ! 00477 ShapeVolume=pi*shape%size_param(1)**2 00478 ELSE 00479 ShapeVolume=4d0*product(shape%size_param(shape%xy(1:2))) 00480 END IF 00481 ELSE 00482 ShapeVolume=2d0*pi*shape%size_param(1)**2*shape%size_param(3) 00483 END IF 00484 CASE (RECTANGULAR_PRISM) 00485 ShapeVolume=2d0**nDim*product(shape%size_param(shape%xy(1:nDim))) 00486 CASE (ELLIPSOID) 00487 IF (nDim == 2) THEN 00488 ShapeVolume=pi*product(shape%size_param(shape%xy(1:2))) 00489 ELSE 00490 ShapeVolume=4d0/3d0*pi*product(shape%size_param(:)) 00491 END IF 00492 CASE (ELLIPTICAL_PRISM) 00493 IF (nDim == 2) THEN 00494 IF (shape%xy(3) == 3) THEN ! 00495 ShapeVolume=pi*product(shape%size_param(shape%xy(1:2))) 00496 ELSE 00497 ShapeVolume=4d0*product(shape%size_param(shape%xy(1:2))) 00498 END IF 00499 ELSE 00500 ShapeVolume=2d0*pi*product(shape%size_param(1:3)) 00501 END IF 00502 END SELECT 00503 END FUNCTION ShapeVolume 00504 00505 END MODULE Shapes 00506