Scrambler  1
FieldLoopRestart/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 FieldLoopRestart 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 
00036 MODULE Problem
00037 
00038   USE DataDeclarations
00039   USE GlobalDeclarations
00040   USE PhysicsDeclarations
00041   USE Ambients
00042   IMPLICIT NONE
00043   SAVE
00044   PRIVATE
00045   
00046   PUBLIC ProblemModuleInit, ProblemGridInit, &
00047          ProblemBeforeStep, ProblemAfterStep, ProblemSetErrFlag, ProblemBeforeGlobalStep
00048 
00049   REAL(KIND=qPrec) :: rho,p,v(3),Ao,thickness,R,phi,theta
00050   REAL(KIND=qPREC), DIMENSION(:,:,:), ALLOCATABLE :: infoq
00051   INTEGER :: sample_res=16
00052 
00053   LOGICAL :: lCooling
00054   LOGICAL :: lResolutionTest=.false.
00055   INTEGER :: iCooling
00056   
00057 
00058 CONTAINS
00059 
00061   SUBROUTINE ProblemModuleInit
00062     INTEGER :: iErr
00063     TYPE(AmbientDef), POINTER :: Ambient
00064     NAMELIST /ProblemData/ rho,p,v,Ao,R,thickness,phi,theta,lCooling,iCooling, lResolutionTest
00065 
00066     OPEN(UNIT=PROBLEM_DATA_HANDLE, FILE='problem.data')
00067        
00068     READ(PROBLEM_DATA_HANDLE,NML=ProblemData)
00069     
00070     CLOSE(PROBLEM_DATA_HANDLE, IOSTAT=iErr)
00071     
00072     
00073     IF (nDim == 2 .AND. ((phi /= zero) .OR. (theta /= zero))) THEN
00074        PRINT *, "CalcConstantProblem() error: phi and theta need to be 0 in a 2D problem."
00075        STOP
00076     END IF
00077     NULLIFY(Ambient)
00078     CALL CreateAmbient(Ambient)
00079     Ambient%density=rho
00080     Ambient%pressure=p
00081     Ambient%velocity=v
00082 END SUBROUTINE ProblemModuleInit
00083  
00086   SUBROUTINE ProblemGridInit(Info)
00087         !! @brief Initializes the grid data according to the requirements of the problem.
00088         !! @param Info A grid structure.        
00089     TYPE (InfoDef) :: Info
00090     INTEGER :: i,j,k,l,m,ii,jj,kk
00091         INTEGER :: rmbc,zrmbc
00092         INTEGER :: mx, my, mz
00093     REAL(KIND=qPrec) :: pos(3), mini_pos(3), x_rotated(3), s, temp(3)
00094     REAL(KIND=qPrec), DIMENSION(:,:,:,:), ALLOCATABLE :: A
00095         INTEGER :: iErr
00096         
00097     
00098 !    IF (ALLOCATED(A))  DEALLOCATE(A)
00099 
00100         ! Calculating the number of ghost cells on each side of the grid.
00101         rmbc = levels(Info%level)%gmbc(1)!CoarsenRatio(Info%level-1) * mbc
00102     SELECT CASE(nDim)
00103     CASE(2)
00104        zrmbc=0
00105     CASE(3)
00106        zrmbc=rmbc
00107     END SELECT
00108 
00109         mx = Info%mX(1)
00110         my = Info%mX(2)
00111         mz = Info%mX(3)
00112 !       write(*,*) Info%mX, Info%xBounds
00113 !        write(*,*) Info%level
00114 !        write(*,*) levels(Info%level)%dx
00115     CALL ConvertTotalToInternalEnergy(Info%q(1:Info%mX(1), 1:Info%mX(2), 1:Info%mX(3), :))
00116 
00117     IF (nDim == 3) THEN
00118            
00119        ALLOCATE(A(1-rmbc:mx+1+rmbc,1-rmbc:my+1+rmbc, 1-zrmbc:mz+1+zrmbc,3), STAT=iErr)
00120            
00121            IF (iErr /= 0) THEN
00122                PRINT *, "DomainInitProblem() error: u unable to allocate 3D array."
00123                    STOP
00124            END IF
00125            
00126            
00127        A=0.0    ! A = Magnetic potential (B = curl(A)).
00128            
00129        DO k=1-zrmbc, mz+1+zrmbc
00130           DO j=1-rmbc, my+1+rmbc
00131              DO i=1-rmbc, mx+1+rmbc
00132                 DO ii=0,0
00133                    DO jj=0,0
00134                       DO kk=0,0
00135                          IF (maxval(abs(A(i,j,k,:))) == 0) THEN
00136 !                            write(*,*) shape(pos)
00137 !                            write(*,*) shape(Info%xBounds)
00138 !                            write(*,*) pos
00139 !                            write(*,*) Info%xBounds(1:3,1)
00140                             pos=Info%xBounds(1:3,1)
00141                             pos=levels(Info%level)%dx
00142                             pos=(/REAL(i,KIND=qPREC)-1d0, REAL(j,KIND=qPREC)-1d0, REAL(k,KIND=qPREC)-1d0/) * levels(Info%level)%dx
00143                             pos=(/REAL(i,8)-1, REAL(j,8)-1, REAL(k,8)-1/) * levels(Info%level)%dx
00144 
00145                                                  
00146                             pos=Info%xBounds(1:3,1)+(/REAL(i,8)-1, REAL(j,8)-1, REAL(k,8)-1/) * levels(Info%level)%dx
00147                             pos(1)=pos(1)+(GxBounds(1,2)-GxBounds(1,1))*(modulo(ii+1,3)-1)
00148                             pos(2)=pos(2)+(GxBounds(2,2)-GxBounds(2,1))*(modulo(jj+1,3)-1)
00149                             pos(3)=pos(3)+(GxBounds(3,2)-GxBounds(3,1))*(modulo(kk+1,3)-1)
00150                                                         
00151                             DO l=1,3
00152                                mini_pos=pos
00153                                DO m=1,sample_res
00154                                   mini_pos(l) = pos(l) + levels(Info%level)%dx*(REAL(m,8)-half)/REAL(sample_res,8) !sub sample each edge
00155                                   x_rotated=rotate_y(rotate_z(mini_pos,-phi),-theta) !rotate coordinates again so that vector along jet-axis is along x-axis
00156                                   s=sqrt(x_rotated(1)**2+x_rotated(2)**2)
00157                                   IF (abs(x_rotated(3)) < half*thickness .AND. s <= R) THEN !inside of loop region
00158                                      temp=emf_source_3D(s)
00159                                      A(i,j,k,l)=A(i,j,k,l)+temp(l)
00160                                   END IF
00161                                END DO
00162                             END DO
00163                          END IF
00164                       END DO
00165                    END DO
00166                 END DO
00167              END DO
00168           END DO
00169        END DO
00170            
00171        A=A/REAL(sample_res, 8)
00172 
00173            ! Calculate magnetic fluxes from A (B = curl(A)).
00174 
00175 IF (MaintainAuxArrays) THEN
00176        Info%aux=0.0
00177 
00178        Info%aux(1-rmbc:mx+rmbc+1, 1-rmbc:my+rmbc, 1-zrmbc:mz+zrmbc, 1)=&
00179              (A(1-rmbc:mx+rmbc+1, 2-rmbc:my+rmbc+1,1-zrmbc:mz+zrmbc, 3) &
00180             - A(1-rmbc:mx+rmbc+1, 1-rmbc:my+rmbc,1-zrmbc:mz+zrmbc, 3))/levels(Info%level)%dX &
00181             - (A(1-rmbc:mx+rmbc+1, 1-rmbc:my+rmbc,2-zrmbc:mz+zrmbc+1, 2)&
00182             - A(1-rmbc:mx+rmbc+1, 1-rmbc:my+rmbc,1-zrmbc:mz+zrmbc, 2))/levels(Info%level)%dx
00183 
00184        Info%aux(1-rmbc:mx+rmbc, 1-rmbc:my+rmbc+1, 1-zrmbc:mz+zrmbc, 2)=&
00185               (A(1-rmbc:mx+rmbc, 1-rmbc:my+rmbc+1,2-zrmbc:mz+zrmbc+1, 1)&
00186             - A(1-rmbc:mx+rmbc, 1-rmbc:my+rmbc+1,1-zrmbc:mz+zrmbc, 1))/levels(Info%level)%dx&
00187             - (A(2-rmbc:mx+rmbc+1, 1-rmbc:my+rmbc+1,1-zrmbc:mz+zrmbc, 3)&
00188             - A(1-rmbc:mx+rmbc, 1-rmbc:my+rmbc+1,1-zrmbc:mz+zrmbc, 3))/levels(Info%level)%dx
00189 
00190        Info%aux(1-rmbc:mx+rmbc, 1-rmbc:my+rmbc, 1-zrmbc:mz+zrmbc+1, 3)=&
00191             (A(2-rmbc:mx+rmbc+1, 1-rmbc:my+rmbc,1-zrmbc:mz+zrmbc+1, 2)&
00192             - A(1-rmbc:mx+rmbc, 1-rmbc:my+rmbc,1-zrmbc:mz+zrmbc+1, 2))/levels(Info%level)%dx&
00193             - (A(1-rmbc:mx+rmbc, 2-rmbc:my+rmbc+1,1-zrmbc:mz+zrmbc+1, 1)&
00194             - A(1-rmbc:mx+rmbc, 1-rmbc:my+rmbc,1-zrmbc:mz+zrmbc+1, 1))/levels(Info%level)%dx
00195 END IF
00196 !       Info%q=0d0
00197 
00198            ! Calculate B-fields
00199 IF (MaintainAuxArrays) THEN
00200        forall(i=1-rmbc:mx+rmbc,j=1-rmbc:my+rmbc, k=1-zrmbc:mz+zrmbc)
00201           Info%q(i,j,k,iBx)=half*(Info%aux(i,j,k,1)+Info%aux(i+1,j,k,1))
00202           Info%q(i,j,k,iBy)=half*(Info%aux(i,j,k,2)+Info%aux(i,j+1,k,2))
00203           Info%q(i,j,k,iBz)=half*(Info%aux(i,j,k,3)+Info%aux(i,j,k+1,3))
00204        end forall
00205 END IF
00206 
00207 !       v = rotate_z(rotate_y(v,theta),phi)
00208     ELSE
00209         
00210                 ! 2D case.
00211 
00212        ALLOCATE(A(1-rmbc:mx+rmbc+1,1-rmbc:my+rmbc+1, 1-zrmbc:mz+zrmbc,1), STAT=iErr)
00213            
00214            IF (iErr /= 0) THEN
00215             PRINT *, "DomainInitProblem() error: unable to allocate 2D array."
00216                    STOP
00217            END IF
00218            
00219        A=0.0            ! A = magnetic potential.
00220            
00221            ! Initialize A.
00222        DO j=1-rmbc, my+rmbc+1
00223           DO i=1-rmbc, mx+rmbc+1
00224              DO ii=0,2 
00225                 DO jj=0,2
00226                    IF (maxval(abs(A(i,j,1,:))) == 0) THEN
00227                       pos=Info%xBounds(1:3,1)+(/REAL(i,8)-1, REAL(j,8)-1, 0d0/) * levels(Info%level)%dx
00228                                           pos(1)=pos(1)+(GxBounds(1,2)-GxBounds(1,1))*(modulo(ii+1,3)-1)
00229                                           pos(2)=pos(2)+(GxBounds(2,2)-GxBounds(2,1))*(modulo(jj+1,3)-1)
00230                       mini_pos=pos
00231                       x_rotated=mini_pos
00232                       s=sqrt(x_rotated(1)**2+x_rotated(2)**2)
00233                       IF (s <= R) THEN !inside of loop region
00234                          temp=emf_source_3D(s)
00235                          A(i,j,1,1)=temp(3)
00236                       END IF
00237                    END IF
00238                 END DO
00239              END DO
00240           END DO
00241        END DO
00242 
00243                 ! Calculate magnetic fluxes (B = curl(A)).
00244 IF (MaintainAuxArrays) THEN
00245        Info%aux=0.0
00246 
00247        Info%aux(1-rmbc:mx+1+rmbc, 1-rmbc:my+rmbc, 1-zrmbc:mz+zrmbc, 1)=&
00248              (A(1-rmbc:mx+rmbc+1, 2-rmbc:my+rmbc+1,1-zrmbc:mz+zrmbc, 1) &
00249             -A(1-rmbc:mx+rmbc+1, 1-rmbc:my+rmbc,1-zrmbc:mz+zrmbc, 1))/levels(Info%level)%dx
00250 
00251        Info%aux(1-rmbc:mx+rmbc, 1-rmbc:my+rmbc+1, 1-zrmbc:mz+zrmbc, 2)=&
00252             -(A(2-rmbc:mx+rmbc+1, 1-rmbc:my+rmbc+1,1-zrmbc:mz+zrmbc, 1)&
00253             -A(1-rmbc:mx+rmbc, 1-rmbc:my+rmbc+1,1-zrmbc:mz+zrmbc, 1))/levels(Info%level)%dx
00254 END IF
00255 
00256 !       Info%q=0d0
00257 
00258        ! Calculate magnetic field values from fluxes.
00259 IF (MaintainAuxArrays) THEN
00260        forall(i=1-rmbc:mx+rmbc,j=1-rmbc:my+rmbc, k=1-zrmbc:mz+zrmbc)
00261           Info%q(i,j,k,iBx)=half*(Info%aux(i,j,k,1)+Info%aux(i+1,j,k,1))
00262           Info%q(i,j,k,iBy)=half*(Info%aux(i,j,k,2)+Info%aux(i,j+1,k,2))
00263        end forall
00264 END IF
00265 
00266     END IF
00267 
00268     ! Set other variables.
00269 
00270 !    Info%q(:,:,:,1)=rho
00271 !    Info%q(:,:,:,2)=rho*v(1)
00272 !    Info%q(:,:,:,3)=rho*v(2)
00273 !    Info%q(:,:,:,4)=rho*v(3)
00274         
00275         ! Calculate the total energy (E_total = E_thermal + E_kinetic + E_magnetic).
00276 
00277     CALL ConvertInternalToTotalEnergy(Info%q(1:Info%mX(1), 1:Info%mX(2), 1:Info%mX(3), :))
00278 
00279     DEALLOCATE(A)
00280 
00281  END SUBROUTINE ProblemGridInit
00282 
00285  SUBROUTINE ProblemBeforeStep(Info)
00286     !! @brief Performs any tasks required before the advance step.
00287     !! @param Info A grid structure.    
00288     TYPE (InfoDef) :: Info
00289     INTEGER, DIMENSION(0:2) :: Steps = (/1,4,11/)
00290     LOGICAL, DIMENSION(0:2) :: RestartTriggered=(/.false.,.false.,.false./)
00291     INTEGER :: i
00292     IF (Info%level > 2) RETURN
00293     i = levels(Info%level)%CurrentLevelStep
00294 
00295     IF (steps(Info%level) == i .AND. MPI_ID == 0 .AND. .NOT. RestartTriggered(Info%level)) THEN
00296        write(*,*) 'Processor ', MPI_ID, ' purposely requesting restart on level ', Info%level, 'step ', i, 'to test code'
00297        lRequestRestart=.true.
00298        RestartTriggered(Info%level)=.true.
00299     END IF
00300  END SUBROUTINE ProblemBeforeStep
00301 
00304   SUBROUTINE ProblemAfterStep(Info)
00305      !! @brief Performs any post-step corrections that are required.
00306      !! @param Info A grid structure.   
00307      TYPE (InfoDef) :: Info
00308      REAL(KIND=qPREC) :: dx, t, x, y, s, Bx, By, err
00309      INTEGER :: i,j
00310      IF (lResolutionTest) THEN
00311         t=levels(Info%level)%tnow+levels(Info%level)%dt
00312         IF (abs(t-.1d0) < 1e-6) THEN
00313            write(*,*) 'saving at t=', t
00314            ALLOCATE(infoq(1:Info%mX(1),1:Info%mX(2),2))
00315            infoq=Info%q(1:Info%mX(1),1:Info%mX(2),1,iBx:iBy)
00316         ELSE IF (abs(t-1.1) < 1e-6) THEN
00317            write(*,*) 'checking at t=', t
00318            dx=levels(Info%level)%dx
00319            err=0
00320            DO i=1, Info%mX(1)
00321               x = GxBounds(1,1)+modulo(Info%xBounds(1,1)+(REAL(i,8)+half)*dx - v(1)*(t-.1) - GxBounds(1,1), GxBounds(1,2)-GxBounds(1,1))
00322               DO j=1, Info%mX(2)
00323                  y = GxBounds(2,1)+modulo(Info%xBounds(2,1)+(REAL(j,8)+half)*dx - v(2)*(t-.1) - GxBounds(2,1), GxBounds(2,2)-GxBounds(2,1))
00324                  s=sqrt(x**2+y**2)
00325 !                 IF (s <= R) THEN 
00326 !                    Bx=-Ao*y/s
00327 !                    By=Ao*x/s
00328 !                 ELSE
00329 !                    Bx=0
00330 !                    By=0
00331 !                 END IF
00332                  if (r > R/2d0) THEN
00333                     err=err+sum(abs(infoq(i,j,:)-info%q(i,j,1,iBx:iBy)))*dx**2
00334                  END if
00335               END DO
00336            END DO
00337            write(*,*) 'resolution, error= ', dx, err
00338         END IF
00339      END IF
00340   END SUBROUTINE ProblemAfterStep
00341 
00344   SUBROUTINE ProblemSetErrFlag(Info)
00345         !! @brief Sets error flags according to problem-specific conditions..
00346         !! @param Info A grid structure.        
00347     TYPE (InfoDef) :: Info
00348   END SUBROUTINE ProblemSetErrFlag
00349 
00350   SUBROUTINE ProblemBeforeGlobalStep(n)
00351      INTEGER :: n
00352   END SUBROUTINE ProblemBeforeGlobalStep
00353 
00354 
00357   FUNCTION emf_source_3D(s)
00358     !! @param s A double-precision number.
00359     REAL(KIND=qPrec) :: s,w(3),emf_source_3D(3)
00360     w=(/0d0,0d0,1d0/)*Ao*(R-s)    
00361     emf_source_3D(1:3)=rotate_z(rotate_y(w,theta),phi) !rotate coordinates again so that vector along jet-axis is along x-axis
00362   END FUNCTION emf_source_3D
00363 
00364 
00365 END MODULE Problem
 All Classes Files Functions Variables