Scrambler  1
cameras.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 !    cameras.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 !#########################################################################
00023 MODULE Cameras
00024    USE GlobalDeclarations
00025    USE PhysicsDeclarations
00026    USE CommonFunctions
00027    IMPLICIT NONE
00028 
00029    TYPE CameraDef
00030       REAL(KIND=qPREC), DIMENSION(3) :: pos
00031       REAL(KIND=qPREC), DIMENSION(3) :: UpVector
00032       REAL(KIND=qPREC), DIMENSION(3) :: Focus
00033       REAL(KIND=qPREC) :: FOV(2)
00034       REAL(KIND=qPREC), DIMENSION(3,3) :: matrix
00035       INTEGER :: res
00036       INTEGER :: iD
00037    END type CameraDef
00038    INTEGER, SAVE :: nCameras=0
00039 
00040 
00041 CONTAINS
00042 
00043    SUBROUTINE CreateCamera(Camera)
00044       TYPE(CameraDef), POINTER :: Camera
00045       ALLOCATE(Camera)
00046       Camera%Focus=half*(SUM(GxBounds,2))
00047       Camera%FOV=(/30d0*Pi/180d0,30d0*Pi/180d0/)
00048       Camera%UpVector=(/0d0,0d0,1d0/)
00049       Camera%pos=Camera%Focus
00050       Camera%pos(2)=Camera%Focus(2)-1d0*(maxval((GxBounds((/1,3/),2)-GxBounds((/1,3/),1))/tan(Camera%FOV)))
00051       CALL UpdateMatrix(Camera)
00052       Camera%res=levels(MaxLevel)%mX(1)
00053       nCameras=nCameras+1
00054       Camera%iD=nCameras
00055    END SUBROUTINE CreateCamera
00056 
00057 
00058    SUBROUTINE UpdateMatrix(Camera)
00059       TYPE(CameraDef) :: Camera
00060       REAL(KIND=qPREC), DIMENSION(3) :: camera_normal(3), ex(3), ey(3)
00061       camera_normal=Camera%focus-Camera%pos
00062       camera_normal=camera_normal/sqrt(sum(camera_normal**2))
00063       ey=Camera%UpVector-SUM(Camera_normal*Camera%UpVector)*Camera_normal
00064       ey=ey/sqrt(sum(ey**2))
00065       ex=Cross3D(camera_normal, ey)
00066       Camera%matrix(:,1)=ex
00067       Camera%matrix(:,2)=ey
00068       Camera%matrix(:,3)=camera_normal
00069 !      write(*,*) 'camera_nomral = ', camera_normal
00070 !      write(*,*) 'ey = ', ey
00071 !      write(*,*) 'ex = ', ex
00072 !      write(*,'(3E25.16)') camera%matrix
00073    END SUBROUTINE UpdateMatrix
00074 
00076 
00077    FUNCTION GetPixel(Camera, pos)
00078       TYPE(CameraDef), POINTER :: Camera
00079       REAL(KIND=qPREC) :: pos(3), GetPixel(2), vec(3)
00080       INTEGER :: i
00081       DO i=1,3
00082          vec(i)=SUM(Camera%matrix(:,i)*(pos-Camera%pos))
00083       END DO
00084 !      write(*,*) 'vec=', vec
00085       GetPixel(1)=atan(vec(1)/vec(3))/Camera%FOV(1)+half
00086       GetPixel(2)=atan(vec(2)/vec(3))/Camera%FOV(2)+half
00087 !      write(*,*) 'atan=', atan(vec(1)/vec(3)), atan(vec(2)/vec(3))
00088    END FUNCTION GetPixel
00089 
00090 
00092    FUNCTION GetRay(Camera, pixel)
00093       TYPE(CameraDef), POINTER :: Camera
00094       REAL(KIND=qPREC) :: GetRay(3), camera_vec(3), pixel(2)
00095       INTEGER :: i
00096       camera_vec(3)=1d0
00097       camera_vec(1)=tan((pixel(1)-half)*Camera%FOV(1))*camera_vec(3)
00098       camera_vec(2)=tan((pixel(2)-half)*Camera%FOV(2))*camera_vec(3)
00099       
00100       DO i=1,3
00101          GetRay(i)=SUM(Camera%matrix(i,:)*camera_vec(:))
00102       END DO
00103       GetRay=GetRay/sqrt(sum(GetRay**2))
00104 
00105 !      write(*,*) GetPixel(camera, Camera%pos+GetRay), pixel
00106 !      write(*,*) Camera%matrix
00107 !      STOP
00108 
00109    END FUNCTION GetRay
00110 
00111 
00112    SUBROUTINE BinCell(Camera, pos, dx, data, rho)
00113       TYPE(CameraDef), POINTER :: Camera
00114       REAL(KIND=qPREC) :: pos(3),xpos(3)
00115       REAL(KIND=qPREC) :: dx, ddx, xbounds(3,2), my_pixel(2), pixel(2)
00116       REAL(KIND=qPREC), DIMENSION(:,:) :: data
00117       REAL(KIND=qPREC) :: rho,a,intersection(6,3), ray(3), max_distance
00118       INTEGER :: ipos(2), sample_res, i, j, k, npoints, min_pixels(2), max_pixels(2), dim, odim(2), edge
00119       xbounds(:,1)=(pos-half*dx)
00120       xbounds(:,2)=(pos+half*dx)
00121       min_pixels=huge(min_pixels(1))
00122       max_pixels=0
00123       DO i=1,2
00124          DO j=1,2
00125             DO k=1,2
00126                my_pixel=GetPixel(Camera, (/xbounds(1,i), xbounds(2,j), xbounds(3,k)/))*REAL(shape(data), KIND=qPREC)+half
00127                min_pixels=max(1,min(min_pixels, floor(my_pixel)))
00128                max_pixels=min(shape(data),max(max_pixels, ceiling(my_pixel)))
00129             END DO
00130          END DO
00131       END DO
00132       DO i=min_pixels(1), max_pixels(1)
00133          DO j=min_pixels(2), max_pixels(2)
00134             pixel=(REAL((/i,j/),KIND=qPREC)-half)/REAL(shape(data), KIND=qPREC)
00135             Ray=GetRay(Camera, pixel)
00136             npoints=0
00137             DO dim=1,3
00138                DO edge=1,2
00139                   ! Camera%pos(dim)+a*ray(dim)=xbounds(dim,edge)
00140                   a=(xbounds(dim,edge)-Camera%pos(dim))/ray(dim)
00141                   xpos=Camera%pos+a*ray
00142                   odim=modulo((/dim,dim+1/),3)+1
00143                   IF (ALL(xpos(odim) >= xbounds(odim,1) .AND. xpos(odim) < xbounds(odim,2))) THEN
00144                      npoints=npoints+1
00145                      intersection(npoints,:)=xpos
00146                   END IF
00147                END DO
00148             END DO
00149             IF (npoints == 0) CYCLE
00150             max_distance=0d0
00151             DO k=2,npoints
00152                max_distance=max(max_distance, sqrt(sum((intersection(k,:)-intersection(1,:))**2)))
00153             END DO
00154             data(i,j)=data(i,j)+rho*max_distance
00155          END DO
00156       END DO
00157    END SUBROUTINE BinCell
00158 
00159 
00160 ! same thing but we want to bin cells distributed over the softening radius with the exponential tail...
00161 
00162    SUBROUTINE BinParticle(Camera, Particle, data, rho)
00163       USE ParticleDeclarations
00164       TYPE(ParticleDef), POINTER :: Particle
00165       TYPE(CameraDef), POINTER :: Camera
00166       REAL(KIND=qPREC) :: pos(3),xpos(3)
00167       REAL(KIND=qPREC) :: dx, ddx, xlower(3)
00168       REAL(KIND=qPREC), DIMENSION(:,:) :: data
00169       REAL(KIND=qPREC) :: rho, wsum, r2, volume
00170       REAL(KIND=qPREC), ALLOCATABLE, DIMENSION(:,:,:) :: w
00171       INTEGER :: ipos(2), sample_res, i, j, k
00172       dx=8d0*sink_dx !diameter of sphere
00173       xlower=Particle%xloc-half*dx
00174       ! angular resolution of image = FOV/res
00175       ! angular size of object = dx/distance
00176       ! If object is much smaller than 1 pixel - then we should subsample...
00177       ! sample_res = angular resolution / 
00178       ! sample res = 
00179       sample_res=8d0*ceiling( (dx/sqrt(sum((Particle%xloc-Camera%pos)**2)) / (camera%FOV(1)/camera%res)))
00180       write(*,*) sample_res
00181       ddx=dx/REAL(sample_res)
00182       xlower=xlower-half*ddx
00183       ALLOCATE(w(sample_res,sample_res,sample_res))
00184       DO i=1,sample_res
00185          xpos(1)=xlower(1)+ddx*i
00186          DO j=1,sample_res
00187             xpos(2)=xlower(2)+ddx*j
00188             DO k=1,sample_res
00189                xpos(3)=xlower(3)+ddx*k
00190                r2=sum((xpos(1:nDim)-Particle%xloc(1:nDim))**2)
00191                IF (r2 < (half*dx)**2) THEN
00192                   w(i,j,k)=1d0!exp(-r2/(Particle%radius*sink_dx)**2)
00193                ELSE
00194                   w(i,j,k)=0d0
00195                END IF
00196             END DO
00197          END DO
00198       END DO
00199       wsum=sum(w)
00200       volume=4d0/3d0*Pi*(dx/2d0)**3
00201       write(*,*) wsum, volume, sample_res
00202       DO i=1,sample_res
00203          xpos(1)=xlower(1)+ddx*i
00204          DO j=1,sample_res
00205             xpos(2)=xlower(2)+ddx*j
00206             DO k=1,sample_res
00207                xpos(3)=xlower(3)+ddx*k
00208  !              r2=sum((xpos(1:nDim)-Particle%xloc(1:nDim))**2)
00209                CALL BinCell(Camera, xpos, ddx, data, rho/volume*w(i,j,k)/wsum*1e6)
00210             END DO
00211          END DO
00212       END DO
00213       DEALLOCATE(w)
00214    END SUBROUTINE BinParticle
00215 
00216 
00217 
00218    SUBROUTINE BinCell_old(Camera, pos, dx, data, rho)
00219       TYPE(CameraDef), POINTER :: Camera
00220       REAL(KIND=qPREC) :: pos(3),xpos(3)
00221       REAL(KIND=qPREC) :: dx, ddx, xlower(3)
00222       REAL(KIND=qPREC), DIMENSION(:,:) :: data
00223       REAL(KIND=qPREC) :: rho
00224       INTEGER :: ipos(2), sample_res, i, j, k
00225       xlower=pos-half*dx
00226       sample_res=2d0*ceiling(camera%FOV(1)/camera%res / (dx/sqrt(sum(pos-Camera%pos)**2)))
00227 !     write(*,*) sample_res
00228       ddx=dx/REAL(sample_res)
00229       xlower=xlower-half*ddx
00230       DO i=1,sample_res
00231          xpos(1)=xlower(1)+ddx*i
00232          DO j=1,sample_res
00233             xpos(2)=xlower(2)+ddx*j
00234             DO k=1,sample_res
00235                xpos(3)=xlower(3)+ddx*k
00236                ipos(1:2)=nint(GetPixel(Camera, xpos)*shape(data))
00237                IF (ALL(ipos(1:2) >= 1) .AND. ALL(ipos(1:2) <= shape(data))) THEN
00238                   data(ipos(1),ipos(2))=data(ipos(1),ipos(2))+ rho*ddx**3
00239                END IF
00240             END DO
00241          END DO
00242       END DO
00243    END SUBROUTINE BinCell_old
00244 
00245 
00246 
00247    SUBROUTINE InputCamera(filehandle, Camera)
00248       TYPE(CameraDef), POINTER :: Camera
00249       TYPE(CameraDef) :: MyCamera
00250       iNTEGER :: filehandle
00251       NAMELIST /CameraData/ MyCamera
00252       MyCamera=Camera
00253       READ(CAMERA_DATA_HANDLE, NML=CameraData) 
00254       write(*,*) 'read camera id=', MyCamera%iD
00255       Camera=MyCamera
00256       CALL UpdateMatrix(Camera)
00257    END SUBROUTINE InputCamera
00258 
00259 END MODULE Cameras
 All Classes Files Functions Variables