Scrambler
1
|
00001 !######################################################################### 00002 ! 00003 ! Copyright (C) 2003-2012 Department of Physics and Astronomy, 00004 ! University of Rochester, 00005 ! Rochester, NY 00006 ! 00007 ! common_functions.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 CommonFunctions 00033 USE GlobalDeclarations 00034 USE PhysicsDeclarations 00035 IMPLICIT NONE 00036 00037 PUBLIC 00038 00039 INTERFACE AddCurl 00040 MODULE PROCEDURE AddCurl2D, AddCurl3D 00041 END INTERFACE 00042 00043 00044 INTEGER, PARAMETER :: INOSMOOTH = 0, ITANH = 1, IDBL_TANH=2, IDECAYING_EXP=3 00045 INTEGER, PUBLIC, PARAMETER :: NOSOFT = 0, !g ~ r/r^2 for all r 00046 SPLINESOFT = 1, !g ~ r/r^2 for r < r_soft and then goes to 0 as r-> 0 00047 PLUMMERSOFT = 2 !g ~ (r^2+r_soft^2)^(-3/2) r 00048 00049 CONTAINS 00050 00051 00054 00059 SUBROUTINE AddCurl2D(b,a,dx) 00060 REAL(KIND=qPREC), DIMENSION(:,:) :: a 00061 REAL(KIND=qPREC), DIMENSION(:,:,:) :: b 00062 REAL(KIND=qPREC) :: dx 00063 INTEGER :: mx, my 00064 mx=size(a,1)-1 00065 my=size(a,2)-1 00066 b(1:mx+1,1:my,1)=b(1:mx+1,1:my,1)+& 00067 (a(1:mx+1,2:my+1)-a(1:mx+1,1:my))/dx 00068 b(1:mx,1:my+1,2)=b(1:mx,1:my+1,2)-& 00069 (a(2:mx+1,1:my+1)-a(1:mx,1:my+1))/dx 00070 END SUBROUTINE AddCurl2D 00071 00076 SUBROUTINE AddCurl3D(b,a,dx) 00077 REAL(KIND=qPREC), DIMENSION(:,:,:,:) :: a 00078 REAL(KIND=qPREC), DIMENSION(:,:,:,:) :: b 00079 REAL(KIND=qPREC) :: dx 00080 INTEGER :: mx, my, mz 00081 mx=size(a,1)-1 00082 my=size(a,2)-1 00083 mz=size(a,3)-1 00084 b(1:mx+1,1:my,1:mz,1)=b(1:mx+1,1:my,1:mz,1)+& 00085 ((a(1:mx+1,2:my+1,1:mz,3)-a(1:mx+1,1:my,1:mz,3)) & 00086 -(a(1:mx+1,1:my,2:mz+1,2)-a(1:mx+1,1:my,1:mz,2)))/dx 00087 b(1:mx,1:my+1,1:mz,2)=b(1:mx,1:my+1,1:mz,2)+& 00088 ((a(1:mx,1:my+1,2:mz+1,1)-a(1:mx,1:my+1,1:mz,1)) & 00089 -(a(2:mx+1,1:my+1,1:mz,3)-a(1:mx,1:my+1,1:mz,3)))/dx 00090 b(1:mx,1:my,1:mz+1,3)=b(1:mx,1:my,1:mz+1,3)+& 00091 ((a(2:mx+1,1:my,1:mz+1,2)-a(1:mx,1:my,1:mz+1,2)) & 00092 -(a(1:mx,2:my+1,1:mz+1,1)-a(1:mx,1:my,1:mz+1,1)))/dx 00093 END SUBROUTINE AddCurl3D 00094 00095 00097 00098 00101 00102 00103 00104 FUNCTION SmoothFunction(r, ifunc, smooth) 00105 REAL(KIND=qPREC) :: r, smooth, SmoothFunction 00106 INTEGER :: ifunc 00107 SELECT CASE(ifunc) 00108 CASE (INOSMOOTH) 00109 SmoothFunction=1d0 00110 CASE(ITANH) 00111 SmoothFunction=smooth_tanh(r, smooth) 00112 CASE(IDBL_TANH) 00113 SmoothFunction= smooth_dbl_tanh(r, smooth) 00114 CASE(IDECAYING_EXP) 00115 SmoothFunction=exp(-r**2/smooth**2) 00116 CASE DEFAULT 00117 PRINT*, 'Unrecognized smoothing option in common_functions.f90' 00118 STOP 00119 END SELECT 00120 END FUNCTION SmoothFunction 00121 00125 FUNCTION smooth_tanh(r,smooth) 00126 REAL(KIND=qPREC) :: smooth_tanh, r, smooth 00127 IF (r >= 1d0) THEN 00128 smooth_tanh=0 00129 ELSEIF (r < 1d0-smooth) THEN 00130 smooth_tanh=1 00131 ELSE 00132 smooth_tanh=tanh((1d0-r)/smooth)/tanh(1d0) 00133 END IF 00134 END FUNCTION smooth_tanh 00135 00139 FUNCTION smooth_dbl_tanh(r,smooth) 00140 REAL(KIND=qPREC) :: smooth_dbl_tanh, r, smooth 00141 IF (r >= 1) THEN 00142 smooth_dbl_tanh=0d0 00143 ELSEIF (r < 1d0-smooth) THEN 00144 smooth_dbl_tanh=1d0 00145 ELSE 00146 smooth_dbl_tanh=half*(1d0+tanh(((1d0-smooth/2d0)-r)/(smooth/2d0))/tanh(1d0)) 00147 END IF 00148 END FUNCTION smooth_dbl_tanh 00149 00150 00152 00153 00154 00155 FUNCTION SphericalRotation(x, theta, phi) 00156 REAL(KIND=qPREC), DIMENSION(3) :: x, SphericalRotation 00157 REAL(KIND=qPREC) :: theta, phi 00158 SphericalRotation=rotate_z(rotate_x(x,-phi), -theta) 00159 END FUNCTION SphericalRotation 00160 00161 FUNCTION InvSphericalRotation(x, theta, phi) 00162 REAL(KIND=qPREC), DIMENSION(3) :: x, InvSphericalRotation 00163 REAL(KIND=qPREC) :: theta, phi 00164 InvSphericalRotation=rotate_x(rotate_z(x,theta), phi) 00165 END FUNCTION InvSphericalRotation 00166 00170 FUNCTION rotate_x(v,a) 00171 REAL(KIND=qPrec) :: v(3), a, rotate_x(3) 00172 rotate_x = (/ v(1), & 00173 COS(a)*v(2) - SIN(a)*v(3), & 00174 SIN(a)*v(2) + COS(a)*v(3)/) 00175 END FUNCTION rotate_x 00176 00180 FUNCTION rotate_y(v,a) 00181 REAL(KIND=qPrec) :: v(3), a, rotate_y(3) 00182 00183 rotate_y = (/COS(a)*v(1) + SIN(a)*v(3), & 00184 v(2), & 00185 -SIN(a)*v(1) + COS(a)*v(3)/) 00186 END FUNCTION rotate_y 00187 00191 FUNCTION rotate_z(v,a) 00192 REAL(KIND=qPrec) :: v(3), a, rotate_z(3) 00193 rotate_z = (/COS(a)*v(1) - SIN(a)*v(2), & 00194 SIN(a)*v(1) + COS(a)*v(2), & 00195 v(3)/) 00196 END FUNCTION rotate_z 00197 00198 FUNCTION rotate_z2D(v,a) 00199 REAL(KIND=qPrec) :: v(2), a, rotate_z2D(2) 00200 rotate_z2D = (/COS(a)*v(1) - SIN(a)*v(2), & 00201 SIN(a)*v(1) + COS(a)*v(2)/) 00202 END FUNCTION rotate_z2D 00203 00204 00205 !Set up and solve the Least Square Fitting based on observed times 00206 !Times is the input data set, ncase is the number of input cases, nvar is the size of the solution vector. 00207 !sls stores the solution vector in an order of (xy,x,y,1) or (xyz,xy,yz,zx,x,y,z,1) 00208 !errls stores the relative error for the fitting 00209 SUBROUTINE LSSolver(JLs,lsdata,sls,errls) 00210 IMPLICIT NONE 00211 REAL(KIND=qPREC), INTENT(IN), DIMENSION(:) :: lsdata 00212 REAL(KIND=qPREC), INTENT(IN), DIMENSION(:,:) :: JLs 00213 REAL(KIND=qPREC), INTENT(INOUT), DIMENSION(:) :: sls 00214 REAL(KIND=qPREC), OPTIONAL, DIMENSION(:) :: errls 00215 00216 REAL(KIND=qPREC), ALLOCATABLE, DIMENSION(:) :: temp 00217 INTEGER :: i,j,nvar 00218 REAL(KIND=qPREC), DIMENSION(:,:), ALLOCATABLE :: Jlsnew 00219 00220 nvar=size(sls,1) 00221 ALLOCATE(temp(nvar)) 00222 ALLOCATE(Jlsnew(1:nvar,1:nvar)) 00223 00224 ! Jlsnew is transpose(Jls)*Jls 00225 DO j=1,nvar 00226 DO i=1,nvar 00227 Jlsnew(i,j)=dot_product(Jls(:,i),Jls(:,j)) 00228 END DO 00229 END DO 00230 00231 ! Invert Jlsnew 00232 CALL MatrixInvert(Jlsnew,nvar) 00233 00234 ! Find the solution by the explicit equation 00235 DO i=1,nvar 00236 temp(i)=dot_product(Jls(:,i),lsdata(:)) 00237 END DO 00238 DO i=1,nvar 00239 sls(i) = dot_product(Jlsnew(i,:),temp(:)) 00240 END DO 00241 00242 IF (present(errls)) THEN 00243 ! compute the error vector of this predictor 00244 DO i=1,size(errls) 00245 errls(i)=abs(lsdata(i)-dot_product(Jls(i,:),sls(:)))/lsdata(i) 00246 END DO 00247 END IF 00248 DEALLOCATE(Jlsnew,temp) 00249 END SUBROUTINE LSSolver 00250 00251 ! Invert any input matrix A using LU decomposition, return its inversion stored in A. 00252 ! Input n is the size of the matrix A. tiny works as the tolerance, which can 00253 ! be varied according to needs 00254 SUBROUTINE MatrixInvert(A,n) 00255 REAL(KIND=qPREC), INTENT(INOUT), DIMENSION(1:n,1:n) :: A 00256 INTEGER, INTENT(IN) :: n 00257 00258 INTEGER :: i,j,k,imax,ii,ll 00259 INTEGER, DIMENSION(:), ALLOCATABLE :: indx 00260 REAL(KIND=qPREC), PARAMETER :: tiny = 1.0d-20 00261 REAL(KIND=qPREC) :: aamax,summ,dum,d 00262 REAL(KIND=qPREC), DIMENSION(:), ALLOCATABLE :: v 00263 REAL(KIND=qPREC), DIMENSION(:,:), ALLOCATABLE :: temp,temp2 00264 00265 d=1d0; 00266 ALLOCATE(indx(1:n)) 00267 ALLOCATE(v(1:n)) 00268 ALLOCATE(temp(1:n,1:n)) 00269 ALLOCATE(temp2(1:n,1:n)) 00270 00271 temp=0d0; temp2=A 00272 DO i = 1, n 00273 aamax = 0d0 00274 DO j = 1, n 00275 IF(abs(A(i,j)).gt.aamax)THEN 00276 aamax = abs(A(i,j)) 00277 END IF 00278 END DO 00279 IF (aamax.eq.0d0)THEN 00280 PRINT*, 'Found singular matrix in LSMatrixInvert. STOP' 00281 STOP 00282 END IF 00283 v(i) = 1d0/aamax 00284 temp(i,i)=1d0 00285 END DO 00286 00287 DO j = 1, n 00288 DO i = 1, j-1 00289 summ = A(i,j) 00290 summ = summ - dot_product(A(i,1:i-1),A(1:i-1,j)) 00291 A(i,j) = summ 00292 END DO 00293 aamax = 0d0 00294 DO i = j, n 00295 summ = A(i,j) 00296 summ = summ - dot_product(A(i,1:j-1),A(1:j-1,j)) 00297 A(i,j) = summ 00298 dum = v(i)*abs(summ) 00299 IF(dum.ge.aamax)THEN 00300 imax = i 00301 aamax = dum 00302 END IF 00303 END DO 00304 IF(j.ne.imax)THEN 00305 DO k = 1, n 00306 dum = A(imax,k) 00307 A(imax,k) = A(j,k) 00308 A(j,k) = dum 00309 END DO 00310 d = -d 00311 v(imax) = v(j) 00312 END IF 00313 indx(j) = imax 00314 IF(A(j,j).eq.0d0)THEN 00315 A(j,j) = tiny 00316 END IF 00317 IF(j.ne.n)THEN 00318 dum = 1d0/A(j,j) 00319 DO i = j+1, n 00320 A(i,j) = A(i,j)*dum 00321 END DO 00322 END IF 00323 END DO 00324 00325 DO j=1,n 00326 ii = 0 00327 DO i=1,n 00328 ll=indx(i) 00329 summ=temp(ll,j) 00330 temp(ll,j)=temp(i,j) 00331 IF(ii.ne.0)THEN 00332 summ=summ-dot_product(A(i,ii:i-1),temp(ii:i-1,j)) 00333 ELSE IF(summ.ne.0d0)THEN 00334 ii=i 00335 END IF 00336 temp(i,j)=summ 00337 END DO 00338 DO i=n,1,-1 00339 summ=temp(i,j) 00340 summ=summ-dot_product(A(i,i+1:n),temp(i+1:n,j)) 00341 temp(i,j)=summ/A(i,i) 00342 END DO 00343 END DO 00344 00345 A=temp 00346 00347 DEALLOCATE(indx) 00348 DEALLOCATE(v) 00349 DEALLOCATE(temp) 00350 DEALLOCATE(temp2) 00351 00352 RETURN 00353 00354 END SUBROUTINE MatrixInvert 00355 00356 FUNCTION GravityPotential(mass,pos,r_soft,soft_function) 00357 REAL(KIND=qPREC), DIMENSION(3) :: pos 00358 REAL(KIND=qPREC) :: GravityPotential 00359 INTEGER :: soft_function 00360 REAL(KIND=qPREC) :: mass, r_soft, r, r2,theta(3),kmin(3), k(3), fact 00361 INTEGER , PARAMETER :: NMax=32 00362 INTEGER :: i,j,l 00363 00364 ! IF (.NOT. lGravityPeriodic) THEN 00365 r2=sum(pos(1:nDim)**2) 00366 r=sqrt(r2) 00367 SELECT CASE(soft_function) 00368 CASE(NOSOFT) 00369 GravityPotential=-1d0/r 00370 CASE(SPLINESOFT) 00371 GravityPotential=SplinePotential(r,r_soft) 00372 CASE(PLUMMERSOFT) 00373 GravityPotential=PlummerPotential(r,r_soft) 00374 END SELECT 00375 ! ELSE 00376 00377 ! rho = M*(delta_func) 00378 ! rho^(k) = M*1d0 00379 ! phi^(k) = 4*Pi*G*M*1d0/k^2 00380 ! phi = int(phi^{k) e(2Pi ikx) dk) = sum(phi^(k) e^(2Pi ikx)) delta_k**nDim 00381 ! = 4*Pi*G*M*sum(e^{2Pi ikx)/k^2) * delta_k**nDim 00382 ! = 4*Pi*G*M*delta_k**nDim * (sum(cos/k^2) 00383 ! has units of G M l^(2-nDim) 00384 ! so for nDim == 3 has units of GM/l 00385 ! for nDim == 2 has units of GM log(l) 00386 ! for nDim == 1 has units of GM l 00387 00388 ! GravityPotential=0d0 00389 00390 00391 ! kmin(1:nDim)=2d0*Pi/(GxBounds(1:nDim,2)-GxBounds(1:nDim,1)) 00392 ! SELECT CASE(nDim) 00393 ! CASE(1) 00394 ! DO i=1,Nmax 00395 ! k(1)=kmin(1)*i 00396 ! theta(1)=k(1)*pos(1) 00397 ! GravityPotential=GravityPotential+cos(theta(1))/k(1)**2 00398 ! END DO 00399 ! CASE(2) 00400 ! DO i=0,Nmax 00401 ! k(1)=kmin(1)*i 00402 ! theta(1)=k(1)*pos(1) 00403 ! DO j=0,Nmax 00404 ! IF (i == 0 .AND. j == 0) CYCLE 00405 ! IF (i == 0 .OR. j == 0) THEN 00406 ! fact=2d0 00407 ! ELSE 00408 ! fact=4d0 00409 ! END IF 00410 ! k(2)=kmin(2)*j 00411 ! theta(2)=k(2)*pos(2) 00412 ! GravityPotential=GravityPotential+fact*cos(theta(1))*cos(theta(2))/sum(k(1:2)**2) 00413 ! END DO 00414 ! END DO 00415 ! CASE(3) 00416 ! DO i=1,Nmax 00417 ! k(1)=kmin(1)*i 00418 ! theta(1)=k(1)*pos(1) 00419 ! DO j=1,Nmax 00420 ! k(2)=kmin(2)*j 00421 ! theta(2)=k(2)*pos(2) 00422 ! DO l=1,Nmax 00423 ! k(3)=kmin(3)*l 00424 ! theta(3)=k(3)*pos(3) 00425 ! GravityPotential=GravityPotential+cos(sum(theta(1:3)))/sum(k(1:3)**2) 00426 ! END DO 00427 ! END DO 00428 ! END DO 00429 ! END SELECT 00430 ! GravityPotential=GravityPotential*4d0*Pi/product(GxBounds(1:nDim,2)-GxBounds(1:nDim,1)) 00431 ! END IF 00432 GravityPotential=GravityPotential*mass*ScaleGrav 00433 END FUNCTION GravityPotential 00434 00435 00436 FUNCTION GravityForce(mass,pos,r_soft,soft_function) 00437 REAL(KIND=qPREC), DIMENSION(3) :: GravityForce, pos 00438 INTEGER :: soft_function 00439 REAL(KIND=qPREC) :: mass, r_soft, r, r2 00440 SELECT CASE(soft_function) 00441 CASE(NOSOFT) 00442 r2=sum(pos(1:nDim)**2) 00443 r=sqrt(r2) 00444 GravityForce=-pos/(r*r2) 00445 CASE(SPLINESOFT) 00446 ! write(*,*) 'spline r_soft=', r_soft 00447 GravityForce=-SplineSoftening(pos,r_soft) 00448 CASE(PLUMMERSOFT) 00449 ! write(*,*) 'plummer r_soft=', r_soft 00450 GravityForce=-PlummerSoftening(pos,r_soft) 00451 END SELECT 00452 GravityForce=GravityForce*mass*ScaleGrav 00453 END FUNCTION GravityForce 00454 00456 FUNCTION SplineSoftening(pos,r_soft) 00457 REAL(KIND=qPREC) :: r,r_soft,y,pos(3),r2 00458 REAL(KIND=qPREC) :: SplineSoftening(3) 00459 REAL(KIND=qPREC), PARAMETER :: 00460 GSC1=1.06666666666666667d+1, !4*(8/3) 00461 GSC2=-38.4d0, !4*(-48/5) 00462 GSC3=32d0, !4*(8) 00463 GSC4=2.13333333333333333d+1, !4*(16/3) 00464 GSC5=-48d0, !4*(-12) 00465 GSC6=38.4d0, !4*(48/5) 00466 GSC7=-1.06666666666666667d+1, !4*(-8/3) 00467 GSC8=-6.66666666666666667d-2 !4*(-1/60) 00468 INTEGER :: i 00469 IF (GravityDim == 3) THEN 00470 r2=sum(pos(1:nDim)**2) 00471 r=sqrt(r2) 00472 y=r/r_soft 00473 IF (r >= r_soft) THEN 00474 SplineSoftening=pos/(r2*r) 00475 ELSE 00476 y=r/r_soft 00477 IF (r_soft**2*r <= TINY(1d0)) THEN 00478 SplineSoftening=0d0 00479 ELSE 00480 IF (r >= half*r_soft) THEN 00481 SplineSoftening=pos * (y*(GSC4+y*(GSC5+y*(GSC6+y*GSC7)))+GSC8/y**2) / (r_soft**2*r) 00482 ELSE 00483 SplineSoftening=pos * (y*(GSC1+y**2*(GSC2+y*GSC3))) / (r_soft**2*r) 00484 END IF 00485 END IF 00486 END IF 00487 ELSE 00488 SplineSoftening=PlummerSoftening2D(pos,r_soft) 00489 END IF 00490 END FUNCTION SplineSoftening 00491 00493 FUNCTION SplinePotential(r,r_soft) 00494 REAL(KIND=qPREC) :: r,r_soft,y,pos(3),r2 00495 REAL(KIND=qPREC) :: SplinePotential 00496 REAL(KIND=qPREC), PARAMETER :: 00497 GSC0=-14d0/5d0, 00498 GSC1=16d0/3d0, 00499 GSC2=-48d0/5d0, 00500 GSC3=32d0/5d0, 00501 GSC4=32d0/3d0, 00502 GSC5=-16d0, 00503 GSC6=48d0/5d0, 00504 GSC7=-32d0/15d0, 00505 GSC8=1d0/15d0, 00506 GSC9=-16d0/5d0 00507 y=r/r_soft 00508 00509 IF (GravityDim == 3) THEN 00510 IF (r >= r_soft) THEN 00511 SplinePotential=-1d0/r 00512 ELSE 00513 y=r/r_soft 00514 IF (r >= half*r_soft) THEN 00515 SplinePotential=(GSC9+(y**2*(GSC4+y*(GSC5+y*(GSC6+y*GSC7)))+GSC8/y )) / r_soft 00516 ELSE 00517 SplinePotential=(GSC0 + (y**2*(GSC1+y**2*(GSC2+y*GSC3)))) / r_soft 00518 END IF 00519 END IF 00520 ELSE 00521 SplinePotential=PlummerPotential2D(r,r_soft) 00522 END IF 00523 END FUNCTION SplinePotential 00524 00526 FUNCTION PlummerSoftening(pos,r_soft) 00527 REAL(KIND=qPREC) :: r,y,pos(3),r_soft 00528 REAL(KIND=qPREC) :: PlummerSoftening(3) 00529 IF (GravityDim == 3) THEN 00530 PlummerSoftening=pos/((sum(pos**2)+r_soft**2)**(1.5d0)) 00531 ELSE 00532 PlummerSoftening=PlummerSoftening2D(pos,r_soft) 00533 END IF 00534 END FUNCTION PlummerSoftening 00535 00536 00538 FUNCTION PlummerSoftening2D(pos,r_soft) 00539 REAL(KIND=qPREC) :: r,y,pos(3),r_soft 00540 REAL(KIND=qPREC) :: PlummerSoftening2D(3) 00541 PlummerSoftening2D=2d0*pos/(sum(pos(1:nDim)**2)+r_soft**2) 00542 END FUNCTION PlummerSoftening2D 00543 00544 00546 FUNCTION PlummerPotential(r,r_soft) 00547 REAL(KIND=qPREC) :: r,r_soft, PlummerPotential 00548 IF (GravityDim == 3) THEN 00549 PlummerPotential=-1d0/sqrt(r**2+r_soft**2) 00550 ELSE 00551 PlummerPotential=PlummerPotential2D(r,r_soft) 00552 END IF 00553 END FUNCTION PlummerPotential 00554 00555 00557 FUNCTION PlummerPotential2D(r,r_soft) 00558 REAL(KIND=qPREC) :: r,r_soft, PlummerPotential2D 00559 PlummerPotential2D=2d0*log(sqrt(r**2+r_soft**2)/R2DEff)! (*exp(-.5d0))) !)sqrt(r**2+r_soft**2)) 00560 END FUNCTION PlummerPotential2D 00561 00562 00566 FUNCTION Cross2D(a,b) 00567 REAL(KIND=qPREC) :: a(2),b(2),Cross2D 00568 Cross2D=(a(1)*b(2)-a(2)*b(1)) 00569 END FUNCTION Cross2D 00570 00574 FUNCTION Cross3D(a,b) 00575 REAL(KIND=qPREC) :: a(3),b(3),Cross3D(3) 00576 Cross3D=(/a(2)*b(3)-a(3)*b(2),a(3)*b(1)-a(1)*b(3),a(1)*b(2)-a(2)*b(1)/) 00577 END FUNCTION Cross3D 00578 00579 00581 FUNCTION SolveCrossEq(b,c) 00582 REAL(KIND=qPrec), DIMENSION(3) :: b,c,SolveCrossEq 00583 SolveCrossEq = (Cross3D(b,c)+sqrt(sum(b**2)-sum(Cross3D(b,c)**2)/sum(b**2))*b)/sum(b**2) 00584 END FUNCTION SolveCrossEq 00585 00586 00587 FUNCTION fade(r, thickness) 00588 REAL(KIND=qPrec) :: fade, r, thickness 00589 IF (r <= 0d0 .OR. thickness <= 0) THEN 00590 fade = 1d0 00591 ELSE IF (r < 1d0) THEN 00592 fade = tanh((1d0-r)/thickness) 00593 ELSE 00594 fade = 0d0 00595 END IF 00596 ! fade=1d0 00597 END FUNCTION fade 00598 00599 ! Returns the associated Legendre Polynomial at point x with orders l and m 00600 FUNCTION pLgndr(l,m,x) 00601 REAL(KIND=qPrec) :: pLgndr 00602 INTEGER :: l, m, ll 00603 REAL(KIND=qPrec) :: x, pll,pmm,pmmp1,somx2 00604 IF (.NOT.(m>=0.and.m<=l.and.ABS(x)<=1.0))THEN 00605 PRINT*, "wrong input in plgndr" 00606 STOP 00607 END IF 00608 pmm=1.0 00609 IF(m>0) THEN 00610 somx2=SQRT((1.0-x)*(1.0+x)) 00611 pmm=profac(1d0,2d0,m)*somx2**m 00612 END IF 00613 IF(l==m) THEN 00614 plgndr=pmm 00615 ELSE 00616 pmmp1=x*(2*m+1)*pmm 00617 IF(l==m+1) THEN 00618 plgndr=pmmp1 00619 ELSE 00620 DO ll=m+2,l 00621 pll=(x*(2*ll-1)*pmmp1-(ll+m-1)*pmm)/(ll-m) 00622 pmm=pmmp1 00623 pmmp1=pll 00624 END DO 00625 plgndr=pll 00626 END IF 00627 END IF 00628 RETURN 00629 END FUNCTION pLgndr 00630 00631 ! Returns Gamma functions at point xx 00632 FUNCTION GAMMAS(xx) 00633 INTEGER :: i 00634 LOGICAL :: gamma_neg 00635 REAL(KIND=qPrec) :: xx, gammas, tmp, x, summ, temp 00636 REAL(KIND=qPrec) , PARAMETER :: stp=2.5066282746310005, pi=3.141592653589793238 00637 REAL(KIND=qPrec) , PARAMETER, DIMENSION(6) :: coef=(/76.18009172947146, 00638 -86.50532032941677,24.0140982408391,-1.231739572450155, 00639 0.120865097386617e-2,-0.5395239384953e-5/) 00640 IF(xx<0)THEN 00641 IF(xx-REAL(FLOOR(xx))==0.) PRINT*, "infinit gamma function value!" 00642 x=1.-xx 00643 gamma_neg=.true. 00644 ELSE 00645 x=xx 00646 gamma_neg=.false. 00647 END IF 00648 summ=0.0; tmp=x+5.5 00649 tmp=(x+0.5)*LOG(tmp)-tmp 00650 DO i=1,6 00651 summ=summ+coef(i)/(x+REAL(i)) 00652 END DO 00653 temp=EXP(tmp+LOG(stp*(1.000000000190015+summ)/x)) 00654 IF(gamma_neg)THEN 00655 gammas=pi/(Sin(pi*xx)*temp) 00656 ELSE 00657 gammas=temp 00658 END IF 00659 RETURN 00660 END FUNCTION GAMMAS 00661 00662 00663 ! Returns Spherical Bessel j at point x with order n 00664 FUNCTION SBesj(n,x) 00665 INTEGER :: i,ind,n 00666 REAL(KIND=qPrec) , PARAMETER :: TOL=1.0e-12, pi=3.141592653589793238 00667 REAL(KIND=qPrec) :: SBesj,x,k,r 00668 IF(x==0.) x=TOL 00669 SBesj=0. 00670 ind=(n+1)/2;k=sqrt(2./(pi*x)) 00671 DO i =0, ind 00672 SBesj=SBesj+(-1.)**(n-i)*(Sin(x)*(2./x)**(n-2*i)*(GAMMAS(REAL(n-i+1,qPrec))/GAMMAS(REAL(i+1,qPrec)))*BNomial(-0.5-REAL(i,qPrec),REAL(n-2*i,qPrec)) 00673 -Cos(x)*(2./x)**(n+1-2*i)*(GAMMAS(REAL(n-i+1,qPrec))/GAMMAS(REAL(i+1,qPrec)))*i*BNomial(-0.5-REAL(i,qPrec),REAL(n-2*i+1,qPrec))) 00674 END DO 00675 SBesj=SBesj*k*sqrt(pi/(2.*x)) 00676 RETURN 00677 END FUNCTION SBesj 00678 00679 ! Returns Spherical Bessel y at point x with order n 00680 FUNCTION SBesy(n,x) 00681 INTEGER :: n 00682 REAL(KIND=qPrec) :: x,SBesy 00683 SBesy=-SBesj(-n-1,x) 00684 RETURN 00685 END FUNCTION SBesy 00686 00687 ! Returns the Binomial function with any real order 00688 FUNCTION BNomial(n,k) 00689 REAL(KIND=qPrec) :: n,k,BNomial 00690 BNomial=GAMMAS(n+1.)/(GAMMAS(k+1.)*GAMMAS(n-k+1.)) 00691 RETURN 00692 END FUNCTION BNomial 00693 00694 ! Return the product of an arithmetic series. 00695 FUNCTION profac(start,inc,num) 00696 INTEGER :: k,num 00697 REAL(KIND=qPrec) :: profac, start, inc, temp 00698 profac=start; temp=start 00699 IF(num>1)THEN 00700 DO k=1,num-1 00701 temp=temp+inc 00702 profac=profac*temp 00703 END DO 00704 END IF 00705 RETURN 00706 END FUNCTION profac 00707 00709 FUNCTION ConvertCylindrical(pos) 00710 REAL(KIND=qPREC), DIMENSION(3) :: pos, ConvertCylindrical 00711 ConvertCylindrical=(/sqrt(sum(pos(1:2)**2)),GetPhi(pos(1), pos(2)),pos(3)/) 00712 00713 END FUNCTION ConvertCylindrical 00714 00716 FUNCTION ConvertSpherical(pos) 00717 REAL(KIND=qPREC), DIMENSION(3) :: pos, ConvertSpherical 00718 ConvertSpherical(1)=sqrt(sum(pos(1:3)**2)) 00719 ConvertSpherical(2)=acos(pos(3)/ConvertSpherical(1)) 00720 ConvertSpherical(3)=GetPhi(pos(1), pos(2)) 00721 END FUNCTION ConvertSpherical 00722 00724 FUNCTION GetPhi(x, y) 00725 REAL(8) :: GetPhi, x, y 00726 GetPhi=acos(x/sqrt(x**2+y**2)) 00727 IF (y < 0) THEN 00728 GetPhi=2d0*Pi-GetPhi 00729 END IF 00730 END FUNCTION GetPhi 00731 00734 FUNCTION JeansLength(rho, temp) 00735 REAL(KIND=qPREC) :: rho, temp, JeansLength 00736 JeansLength=sqrt(gamma*temp*Pi/ScaleGrav/rho) 00737 END FUNCTION JeansLength 00738 00739 00741 SUBROUTINE random_sphere(A) 00742 REAL(8), DIMENSION(:) :: A 00743 REAL :: rand 00744 REAL(8) :: theta, phi 00745 call random_number(rand) 00746 theta=acos((1d0-2d0*rand)) 00747 call random_number(rand) 00748 phi=rand*2d0*ACOS(-1d0) 00749 A=(/sin(theta)*cos(phi),sin(theta)*sin(phi),cos(theta)/) 00750 END SUBROUTINE random_sphere 00751 00752 SUBROUTINE random_circle(A) 00753 REAL(8), DIMENSION(:) :: A 00754 REAL :: rand 00755 REAL(8) :: phi 00756 call random_number(rand) 00757 phi=rand*2d0*ACOS(-1d0) 00758 A(1:2)=(/cos(phi),sin(phi)/) 00759 END SUBROUTINE random_circle 00760 00761 FUNCTION SpectraK(ipos, mB) 00762 REAL(KIND=qPREC) :: SpectraK(3) 00763 INTEGER :: ipos(:) 00764 INTEGER :: mB(:,:) 00765 INTEGER :: i 00766 SpectraK=0 00767 DO i=1,nDim 00768 IF (ipos(i)==mB(i,1)) THEN 00769 SpectraK(i)=0 00770 ELSEIF (ipos(i) < sum(mB(i,:))/2) THEN 00771 SpectraK(i)=(ipos(i)-mB(i,1))/REAL(mB(i,2)-mB(i,1)+1)*(maxval(mB(:,2)-mB(:,1)+1)) 00772 ELSE 00773 SpectraK(i)=(ipos(i)-mB(i,2)-1)/REAL(mB(i,2)-mB(i,1)+1)*(maxval(mB(:,2)-mB(:,1)+1)) 00774 END IF 00775 END DO 00776 END FUNCTION SpectraK 00777 00778 00779 FUNCTION Curl2D(v, dx) 00780 REAL(KIND=qPREC) :: Curl2D, dx 00781 REAL(KIND=qPREC), DIMENSION(3,3,2) :: v 00782 Curl2D=(v(3,2,2)-v(1,2,2)-v(2,3,1)+v(2,1,1))/(2d0*dx) 00783 END FUNCTION Curl2D 00784 00785 00786 FUNCTION Curl3D(v, dx) 00787 REAL(KIND=qPREC) :: Curl3D(3), dx 00788 REAL(KIND=qPREC), DIMENSION(3,3,3,3) :: v 00789 Curl3D(1)=(v(2,3,2,3)-v(2,1,2,3)-v(2,2,3,2)+v(2,2,1,2))/(2d0*dx) 00790 Curl3D(2)=(v(2,2,3,1)-v(2,2,1,1)-v(3,2,2,3)+v(1,2,2,3))/(2d0*dx) 00791 Curl3D(3)=(v(3,2,2,2)-v(1,2,2,2)-v(2,3,2,1)+v(2,1,2,1))/(2d0*dx) 00792 END FUNCTION Curl3D 00793 00794 00795 00796 END MODULE CommonFunctions