Scrambler  1
outflows.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 !    outflows.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 
00036 
00039 MODULE Outflows
00040    USE GlobalDeclarations
00041    USE DataDeclarations
00042    USE PhysicsDeclarations
00043    USE CommonFunctions
00044    USE DataInfoOps
00045    USE EOS
00046    USE ObjectDeclarations
00047    IMPLICIT NONE
00048    private :: emf_source_3D
00049 
00051    REAL(KIND=qPREC), PARAMETER :: default_buffer=2d0
00052 
00053    TYPE OutflowDef
00054       LOGICAL :: active=.false.
00055       LOGICAL :: initialize=.false.
00056       REAL(KIND=qPREC) :: density=0             !Density of Jet
00057       REAL(KIND=qPREC) :: temperature=0         !Temperature of Jet
00058       REAL(KIND=qPREC) :: velocity=0            !Velocity of Jet
00059       REAL(KIND=qPREC) :: theta=0               !Angle between jet and x-axis 
00060       REAL(KIND=qPREC) :: phi=0                 !3D: Angle around x-axis from xy plane // 2D: angle around z-axis from x-axis
00061       REAL(KIND=qPREC) :: B=0                   !Strength of peak magnetic field in computational units
00062       REAL(KIND=qPREC) :: begin=0               !When does jet turn on
00063       REAL(KIND=qPREC) :: duration=0            !How long does jet run for
00064       REAL(KIND=qPREC) :: decay=0               !How quickly does the jet decay
00065       REAL(KIND=qPREC) :: radius=0              !Radius of jet launching region
00066       REAL(KIND=qPREC) :: thickness=0           !Thickness of source region
00067       REAL(KIND=qPREC) :: open_angle=0          !Jet opening angle (half angle) in radians
00068       REAL(KIND=qPREC) :: oscAmplitude=0        !Amplitude of velocity Oscillations
00069       REAL(KIND=qPREC) :: oscPeriod=0           !Period of velocity oscillations
00070       REAL(KIND=qPREC) :: oscPhase=0            !Initial phase of velocity oscillations
00071       REAL(KIND=qPREC) :: PrecessAmplitude=0    !Precession angle in radians - measured from x' axis
00072       REAL(KIND=qPREC) :: PrecessPeriod=0       !Period of precession
00073       REAL(KIND=qPREC) :: PrecessPhase=0        !Initial phase of precession - measured from x'y' plane in rotated jet frame.
00074       REAL(KIND=qPREC) :: mu=0                  !Inverse beta for jet magnetic field
00075       REAL(KIND=qPREC) :: finish=0      
00076       REAL(KIND=qPREC) :: buffer=default_buffer !Distance from center of jet to source region
00077       REAL(KIND=qPREC) :: offset=0
00078       REAL(KIND=qPREC) :: r_inner=0
00079       REAL(KIND=qPREC) :: r_outer=0
00080       REAL(KIND=qPREC) :: radius_needed
00081       REAL(KIND=qPREC), DIMENSION(3) :: source_vel=(/0,0,0/)  !Outflow velocity (in direction of outflow axis)
00082       REAL(KIND=qPREC), DIMENSION(3) :: position=(/0,0,0/)  !Outflow location
00083       REAL(KIND=qPREC) :: t0=0
00084       REAL(KIND=qPREC), DIMENSION(3,2) :: xBounds !physical extent of outflow box
00085       REAL(KIND=qPREC) :: mass=0                    !Current mass of outflow object
00086       REAL(KIND=qPREC) :: accretionrate=0                !Current accretion rate
00087       REAL(KIND=qPREC) :: massloss=0                !mass lossed to interior cells of grid
00088       INTEGER :: id
00089       INTEGER :: ObjId
00090 
00091 
00092    END TYPE OutflowDef
00093 
00094    !new declaration
00095    TYPE pOutflowDef
00096       TYPE(OutflowDef), POINTER :: ptr
00097    END TYPE pOutflowDef
00098   TYPE(pOutflowDef) :: pOutflow
00099   !
00100 
00101    INTEGER ::  nOutflowObjects = 0
00102    SAVE
00103 CONTAINS
00104 
00105 
00111 
00112 
00113      SUBROUTINE CreateOutflowObject(Outflow, density, temperature, velocity)
00114     !INTEGER :: dummy
00115     TYPE(OutflowDef),POINTER :: Outflow
00116     REAL(KIND=qPREC), OPTIONAL :: density, temperature, velocity
00117     
00118     ALLOCATE(Outflow)
00119     nOutflowObjects=nOutflowObjects+1
00120     Outflow%id=nOutflowObjects
00121 
00122     IF (Present(density)) OutFlow%density=density
00123     IF (Present(temperature)) Outflow%temperature=temperature
00124     IF (Present(velocity)) OutFlow%velocity=velocity
00125     
00126     CALL UpdateOutflow(Outflow)
00127     CALL AddOutflowObjToList(Outflow)
00128   END SUBROUTINE CreateOutflowObject
00129 
00130 
00133    SUBROUTINE UpdateOutflow(Outflow)
00134       TYPE(OutflowDef), POINTER :: Outflow
00135       CALL SetOutflowGeometry(Outflow)
00136       CALL SetOutflowBounds(Outflow)     
00137       Outflow%finish=outflow%begin+outflow%duration+outflow%decay
00138    END SUBROUTINE UpdateOutflow
00139 
00140    SUBROUTINE AddOutflowObjToList(Outflow)
00141       TYPE(OutflowDef), POINTER :: Outflow
00142       TYPE(ObjectDef), POINTER :: Object
00143       Outflow%ObjId = ObjectListAdd(Object,OutflowOBJ)
00144       pOutflow%ptr => Outflow
00145       len = size(transfer(pOutflow, dummy_char))
00146       ALLOCATE(Object%storage(len))
00147       Object%storage = transfer(pOutflow,Object%storage)
00148    END SUBROUTINE AddOutflowObjToList
00149 
00150   SUBROUTINE DestroyOutflowObject(Outflow,id)
00151       TYPE(OutflowDef),POINTER :: Outflow
00152       TYPE(ObjectDef),POINTER :: Object
00153       INTEGER,OPTIONAL :: id
00154 
00155       IF(PRESENT(id)) THEN
00156         Object => ObjectListFind(id)
00157         IF (ASSOCIATED(Object) .AND. Object%type == OUTFLOWOBJ) THEN
00158           pOutflow = transfer(Object%storage,pOutFlow)
00159           DEALLOCATE(pOutflow%ptr)
00160           NULLIFY(pOutflow%ptr)
00161           CALL ObjectListRemove(id)
00162         ENDIF
00163       ELSE
00164         CALL ObjectListRemove(Outflow%ObjId)
00165         DEALLOCATE(Outflow)
00166         NULLIFY(Outflow)
00167       ENDIF
00168   END SUBROUTINE DestroyOutflowObject
00169 
00170   SUBROUTINE OutflowGridInit(Info, Outflow)
00171       TYPE(InfoDef) :: Info
00172       TYPE(OutflowDef), POINTER :: Outflow
00173       CALL OutflowBeforeStep(Info, Outflow)
00174   END SUBROUTINE OutflowGridInit
00175 
00176    SUBROUTINE OutflowBeforeStep(Info, Outflow)
00177       TYPE(InfoDef) :: Info
00178       TYPE(OutflowDef), POINTER :: Outflow
00179        IF (Outflow%begin < levels(Info%level)%tnow .AND. Outflow%finish > levels(Info%level)%tnow) THEN
00180           CALL SourceOutflow(Info, Outflow, levels(Info%level)%tnow, levels(Info%level)%dt, IEVERYWHERE)
00181        END IF
00182    END SUBROUTINE OutflowBeforeStep
00183 
00184   SUBROUTINE OutflowSetErrFlag(Info, Outflow)
00185       TYPE(InfoDef) :: Info
00186       TYPE(OutflowDef), POINTER :: Outflow
00187       INTEGER, DIMENSION(:,:,:), POINTER :: mSs
00188       REAL(KIND=qPREC), DIMENSION(3) :: offset
00189       REAL(KIND=qPREC), POINTER, DIMENSION(:,:) :: offsets
00190       INTEGER :: nOverlaps, nGhost
00191       INTEGER, DIMENSION(3,2) :: mS
00192       INTEGER :: i
00193       nghost=0
00194          CALL CalcPhysicalOverlaps(Info, Outflow%xBounds, mSs, nOverlaps, offsets, IEVERYWHERE, lHydroPeriodic, nGhost)
00195          DO i=1,nOverlaps
00196             mS=mSs(i,:,:)
00197             Info%ErrFlag(ms(1,1):mS(1,2),mS(2,1):mS(2,2),mS(3,1):ms(3,2))=1
00198          END DO
00199          IF (nOverlaps > 0) THEN
00200             DEALLOCATE(mSs)
00201             NULLIFY(mSs)
00202          END IF
00203    END SUBROUTINE OutflowSetErrFlag
00204 
00205    SUBROUTINE OutflowBeforeGlobalStep(n)
00206       INTEGER :: n
00207    END SUBROUTINE OutflowBeforeGlobalStep
00208 
00209 
00210 
00214    SUBROUTINE PlaceOutflow(Info,Outflow,location)
00215       TYPE(InfoDef) :: Info
00216       Type(OutflowDef) :: Outflow
00217       INTEGER :: i,j,k,n,m,ii,jj,kk,location
00218       INTEGER, DIMENSION(3,2) :: mS    
00219       INTEGER, POINTER, DIMENSION(:,:,:) :: mSs
00220       REAL(KIND=qPREC), DIMENSION(3) :: offset
00221       REAL(KIND=qPREC), POINTER, DIMENSION(:,:) :: offsets
00222       REAL(KIND=qPREC), DIMENSION(3) :: xpos, pos
00223       REAL(KIND=qPREC) :: sample_fact, q_fact, dx
00224       INTEGER :: sample_res(3), nOverlaps
00225       REAL(KIND=qPREC), DIMENSION(:), ALLOCATABLE :: q_source, emf_source
00226       REAL(KIND=qPREC), DIMENSION(:,:,:,:), ALLOCATABLE :: emf
00227       dx=levels(Info%level)%dX
00228       xpos=0          
00229       CALL CalcPhysicalOverlaps(Info, Outflow%xBounds, mSs, nOverlaps, offsets, location,lHydroPeriodic)
00230       IF (nOverlaps > 0) THEN
00231          sample_res=1
00232          sample_res(1:nDim)=min(2**(MaxLevel-Info%level),8)
00233          sample_fact=1d0/REAL(sample_res(1),8)
00234          q_fact=sample_fact**nDim
00235          ALLOCATE(q_Source(NrHydroVars))
00236          DO n=1,nOverlaps
00237             mS=mSs(n,:,:)
00238             offset=offsets(n,:)
00239             !Set up emf first - aux fields first
00240             CALL ConvertTotalToInternalEnergy(Info%q(mS(1,1):mS(1,2), mS(2,1):mS(2,2),mS(3,1):mS(3,2),:))           
00241             IF (MaintainAuxArrays .AND. Outflow%B > 0d0) THEN
00242                IF (nDim == 2) THEN
00243                   ALLOCATE(emf(mS(1,1):mS(1,2)+1, mS(2,1):mS(2,2)+1,1,1))
00244                   emf=0
00245                   DO i=mS(1,1),mS(1,2)+1
00246                      xpos(1)=Info%xBounds(1,1)+offset(1)+(i-1)*dx
00247                      DO j=mS(2,1),mS(2,2)+1
00248                         xpos(2)=Info%xBounds(2,1)+offset(2)+(j-1)*dx
00249                         emf(i,j,1,1)=emf_Launch_2D(Outflow,xpos)                   
00250                      END DO
00251                   END DO
00252                   CALL AddCurl(Info%aux(mS(1,1):mS(1,2)+1,mS(2,1):mS(2,2)+1,1,1:2),emf(:,:,1,1),dx)
00253                   DEALLOCATE(emf)
00254                ELSE IF (nDim == 3) THEN
00255                   ALLOCATE(emf(mS(1,1):mS(1,2)+1, mS(2,1):mS(2,2)+1,mS(3,1):mS(3,2)+1,3), emf_source(3))
00256                   emf=0
00257                   DO i=mS(1,1),mS(1,2)
00258                      xpos(1)=Info%xBounds(1,1)+offset(1)+(i-1)*dx
00259                      DO j=mS(2,1),mS(2,2)
00260                         xpos(2)=Info%xBounds(2,1)+offset(2)+(j-1)*dx
00261                         DO k=mS(3,1),mS(3,2)
00262                            xpos(3)=Info%xBounds(3,1)+offset(3)+(k-1)*dx
00263                            DO m=1,3
00264                               pos=xpos
00265                               DO ii=1,sample_res(1)
00266                                  pos(m)=xpos(m)+(REAL(ii, 8)-half)*dx*sample_fact
00267                                  emf_source=emf_launch_3D(Outflow,pos)
00268                                  emf(i,j,k,m)=emf(i,j,k,m)+emf_source(m)
00269                               END DO
00270                               emf(i,j,k,m)=emf(i,j,k,m)*sample_fact
00271                            END DO
00272                         END DO
00273                      END DO
00274                   END DO
00275                   CALL AddCurl(Info%aux(mS(1,1):mS(1,2)+1,mS(2,1):mS(2,2)+1,mS(3,1):mS(3,2)+1,1:3),emf,dx)
00276                   DEALLOCATE(emf, emf_source)                 
00277                END IF
00278                CALL UpdateAux(Info, mS)          
00279             END IF
00280             ! Now set up cell centered quantities (density and momentum)
00281             DO k=mS(3,1),mS(3,2)
00282                xpos(3)=Info%xBounds(3,1)+offset(3)+REAL(k-1,8)*dx
00283                DO j=mS(2,1),mS(2,2)
00284                   xpos(2)=Info%xBounds(2,1)+offset(2)+REAL(j-1,8)*dx
00285                   DO i=mS(1,1),mS(1,2)
00286                      xpos(1)=Info%xBounds(1,1)+offset(1)+REAL(i-1,8)*dx
00287                      q_Source=0
00288                      DO kk=1,sample_res(3)
00289                         pos(3)=xpos(3)+(REAL(kk, 8)-half)*dx*sample_fact
00290                         DO jj=1,sample_res(2)
00291                            pos(2)=xpos(2)+(REAL(jj, 8)-half)*dx*sample_fact
00292                            DO ii=1,sample_res(1)
00293                               pos(1)=xpos(1)+(REAL(ii, 8)-half)*dx*sample_fact
00294                               CALL OutflowLaunch(Outflow,pos,Info%q(i,j,k,1:NrHydroVars),q_source)
00295                            END DO
00296                         END DO
00297                      END DO
00298                      Info%q(i,j,k,1:NrHydroVars)=Info%q(i,j,k,1:NrHydroVars)+q_source*q_fact                 
00299                      IF (i >= 1 .AND. i <= Info%mX(1) .AND. j >= 1 .AND. j <= Info%mX(2) .AND. k >= 1 .AND. k <= Info%mX(3)) Outflow%massloss=Outflow%massloss+q_source(1)*q_fact
00300                   END DO
00301                END DO
00302             END DO
00303             CALL ConvertInternalToTotalEnergy(Info%q(mS(1,1):mS(1,2), mS(2,1):mS(2,2),mS(3,1):mS(3,2),:))           
00304          END DO
00305          DEALLOCATE(q_Source)
00306       if (associated(mSs)) then
00307         DEALLOCATE(mSs)
00308         nullify(mSs)
00309       endif
00310       if (associated(offsets)) then
00311         DEALLOCATE(offsets)
00312         nullify(offsets)
00313       endif
00314  END IF
00315    END SUBROUTINE PlaceOutflow
00316 
00317 
00321    SUBROUTINE SourceOutflow(Info,Outflow,t,dt,location)
00322       TYPE(InfoDef) :: Info
00323       Type(OutflowDef) :: Outflow
00324       INTEGER :: i,j,k,n,m,ii,jj,kk,location
00325       INTEGER, DIMENSION(3,2) :: mS    
00326       INTEGER, POINTER, DIMENSION(:,:,:) :: mSs
00327       REAL(KIND=qPREC), DIMENSION(3) :: offset
00328       REAL(KIND=qPREC), POINTER, DIMENSION(:,:) :: offsets
00329       REAL(KIND=qPREC), DIMENSION(3) :: xpos, pos
00330       REAL(KIND=qPREC) :: sample_fact, q_fact, dx, t, dt
00331       INTEGER :: sample_res(3), nOverlaps
00332       REAL(KIND=qPREC), DIMENSION(:), ALLOCATABLE :: q_source, emf_source
00333       REAL(KIND=qPREC), DIMENSION(:,:,:,:), ALLOCATABLE :: emf
00334       dx=levels(Info%level)%dX
00335       xpos=0          
00336       CALL CalcPhysicalOverlaps(Info, Outflow%xBounds, mSs, nOverlaps, offsets,location,lHydroPeriodic)
00337       IF (nOverlaps > 0) THEN
00338          sample_res=1
00339          sample_res(1:nDim)=min(2**(MaxLevel-Info%level),8)
00340          sample_fact=1d0/REAL(sample_res(1),8)
00341          q_fact=sample_fact**nDim
00342          ALLOCATE(q_Source(NrHydroVars))
00343          DO n=1,nOverlaps
00344             mS=mSs(n,:,:)
00345             offset=offsets(n,:)
00346             !Set up emf first - aux fields first
00347             CALL ConvertTotalToInternalEnergy(Info%q(mS(1,1):mS(1,2), mS(2,1):mS(2,2),mS(3,1):mS(3,2),:))           
00348             IF (MaintainAuxArrays .AND. Outflow%B > 0d0) THEN
00349                IF (nDim == 2) THEN
00350                   ALLOCATE(emf(mS(1,1):mS(1,2)+1, mS(2,1):mS(2,2)+1,1,1))
00351                   emf=0
00352                   DO i=mS(1,1),mS(1,2)+1
00353                      xpos(1)=Info%xBounds(1,1)+offset(1)+(i-1)*dx
00354                      DO j=mS(2,1),mS(2,2)+1
00355                         xpos(2)=Info%xBounds(2,1)+offset(2)+(j-1)*dx
00356                         emf(i,j,1,1)=emf_Source_2D(Outflow,xpos,t,dt)                   
00357                      END DO
00358                   END DO
00359                   CALL AddCurl(Info%aux(mS(1,1):mS(1,2)+1,mS(2,1):mS(2,2)+1,1,1:2),emf(:,:,1,1),dx)
00360                   DEALLOCATE(emf)
00361                ELSE IF (nDim == 3) THEN
00362                   ALLOCATE(emf(mS(1,1):mS(1,2)+1, mS(2,1):mS(2,2)+1,mS(3,1):mS(3,2)+1,3), emf_source(3))
00363                   emf=0
00364                   DO i=mS(1,1),mS(1,2)
00365                      xpos(1)=Info%xBounds(1,1)+offset(1)+(i-1)*dx
00366                      DO j=mS(2,1),mS(2,2)
00367                         xpos(2)=Info%xBounds(2,1)+offset(2)+(j-1)*dx
00368                         DO k=mS(3,1),mS(3,2)
00369                            xpos(3)=Info%xBounds(3,1)+offset(3)+(k-1)*dx
00370                            DO m=1,3
00371                               pos=xpos
00372                               DO ii=1,sample_res(1)
00373                                  pos(m)=xpos(m)+(REAL(ii, 8)-half)*dx*sample_fact
00374                                  emf_source=emf_source_3D(Outflow,pos,t,dt)
00375                                  emf(i,j,k,m)=emf(i,j,k,m)+emf_source(m)
00376                               END DO
00377                               emf(i,j,k,m)=emf(i,j,k,m)*sample_fact
00378                            END DO
00379                         END DO
00380                      END DO
00381                   END DO
00382                   CALL AddCurl(Info%aux(mS(1,1):mS(1,2)+1,mS(2,1):mS(2,2)+1,mS(3,1):mS(3,2)+1,1:3),emf,dx)
00383                   DEALLOCATE(emf, emf_source)                 
00384                END IF
00385                CALL UpdateAux(Info, mS)          
00386             END IF
00387             ! Now set up cell centered quantities (density and momentum)
00388             DO k=mS(3,1),mS(3,2)
00389                xpos(3)=Info%xBounds(3,1)+offset(3)+REAL(k-1,8)*dx
00390                DO j=mS(2,1),mS(2,2)
00391                   xpos(2)=Info%xBounds(2,1)+offset(2)+REAL(j-1,8)*dx
00392                   DO i=mS(1,1),mS(1,2)
00393                      xpos(1)=Info%xBounds(1,1)+offset(1)+REAL(i-1,8)*dx
00394                      q_Source=0
00395                      DO kk=1,sample_res(3)
00396                         pos(3)=xpos(3)+(REAL(kk, 8)-half)*dx*sample_fact
00397                         DO jj=1,sample_res(2)
00398                            pos(2)=xpos(2)+(REAL(jj, 8)-half)*dx*sample_fact
00399                            DO ii=1,sample_res(1)
00400                               pos(1)=xpos(1)+(REAL(ii, 8)-half)*dx*sample_fact
00401                               CALL OutflowSource(Outflow,pos,Info%q(i,j,k,1:NrHydroVars),q_source,t,dt)
00402                            END DO
00403                         END DO
00404                      END DO
00405                      Info%q(i,j,k,1:NrHydroVars)=Info%q(i,j,k,1:NrHydroVars)+q_source*q_fact                 
00406                      IF (i >= 1 .AND. i <= Info%mX(1) .AND. j >= 1 .AND. j <= Info%mX(2) .AND. k >= 1 .AND. k <= Info%mX(3)) Outflow%massloss=Outflow%massloss+q_source(1)*q_fact
00407                   END DO
00408                END DO
00409             END DO
00410             CALL ConvertInternalToTotalEnergy(Info%q(mS(1,1):mS(1,2), mS(2,1):mS(2,2),mS(3,1):mS(3,2),:))           
00411          END DO
00412       if (associated(mSs)) then
00413         DEALLOCATE(mSs)
00414         nullify(mSs)
00415       endif
00416       if (associated(offsets)) then
00417         DEALLOCATE(offsets)
00418         nullify(offsets)
00419       endif
00420    END IF
00421    END SUBROUTINE SourceOutflow
00422 
00423 
00424 
00425    SUBROUTINE OutflowLaunch(Outflow, pos, q, dq)
00426       REAL(KIND=qPrec), DIMENSION(:) :: pos,q,dq
00427       REAL(KIND=qPrec) :: s, r, angle, w(3), vddr, pr, pT, mag_press,therm_press,rz,rz2,s2,x,t,dt,rho_Outflow,f,b_curr(3),precessPhase,oscPhase,velocity,rpos(3)
00428       TYPE(OutflowDef) :: Outflow
00429       pos=pos-Outflow%position !refere to the relative position of source
00430       dq=0
00431       precessphase=0
00432       oscphase=0
00433 
00434       IF (Outflow%precessperiod > 0) PrecessPhase=Outflow%PrecessPhase
00435       IF (Outflow%oscPeriod > 0) oscPhase=Outflow%oscPhase
00436       velocity=Outflow%velocity+sin(oscPhase)*Outflow%oscAmplitude
00437       !    write(*,*) t,dt,precessphase,oscphase,velocity
00438       SELECT CASE(nDim)
00439       CASE(2)
00440 
00441          rpos(1:2)=rotate_z2D(pos(1:2),-Outflow%phi-PrecessPhase) !rotate to coordinates in Outflow frame
00442          x=Outflow%offset+abs(rpos(1)) !distance to 'center' for Outflows
00443          s2=rpos(2)**2
00444          s=abs(rpos(2))
00445 
00446          IF (Outflow%open_angle == 0) THEN
00447             IF (s < Outflow%radius .AND. x < Outflow%r_outer) THEN
00448                IF (Outflow%B > 0) THEN
00449                   write(*,*) "Magnetized Outflows in 2D not supported"
00450                   STOP
00451                END IF
00452                dq(1)=Outflow%density-q(1)
00453                therm_press=Outflow%density*Outflow%temperature
00454                IF (iE .ne. 0) dq(iE)=therm_press-q(iE)
00455                f=fade(s/Outflow%radius, .5d0)
00456                vddr=f*velocity/Outflow%thickness
00457                IF (x > Outflow%buffer) THEN
00458                   dq(ivx:ivy)=Outflow%density*(vddr*(x-Outflow%buffer)*rotate_z2D((/sign(1d0,rpos(1)), 0d0/),Outflow%phi+precessPhase)+Outflow%source_vel(1:2))-q(ivx:ivy) !rotate coordinates again so that vector along Outflow-axis is along x-axis
00459                ELSE
00460                   dq(ivx:ivy)=Outflow%density*(Outflow%source_vel(1:2))-q(ivx:ivy)
00461                END IF
00462             END IF
00463          ELSEIF (abs(Outflow%open_angle-half*Pi) < .01) THEN
00464             !          angle=atan(s/x)
00465             r=sqrt(s2+x**2)
00466             IF (r < Outflow%r_outer) THEN
00467                dq(1)=Outflow%density-q(1)
00468                IF (iE .ne. 0) dq(iE)=Outflow%density*Outflow%temperature-q(iE)
00469                IF (r > Outflow%r_inner) THEN
00470                   vddr=velocity/Outflow%thickness
00471                   dq(ivx:ivy)=Outflow%density*(vddr*(r-Outflow%r_inner)*rotate_z2D(rpos(1:2)/r,Outflow%phi+precessPhase)+Outflow%source_vel(1:2))-q(ivx:ivy) !rotate coordinates again so that vector along Outflow-axis is along x-axis
00472                ELSE
00473                   dq(ivx:ivy)=Outflow%density*(Outflow%source_vel(1:2))-q(ivx:ivy)
00474                END IF
00475             END IF
00476          ELSE
00477             angle=atan(s/x)
00478             r=sqrt(s2+x**2)
00479             IF (angle < Outflow%open_angle .AND. r < Outflow%r_outer) THEN
00480                dq(1)=Outflow%density-q(1)
00481                IF (iE .ne. 0) dq(iE)=Outflow%density*Outflow%temperature-q(iE)
00482                IF (r > Outflow%r_inner) THEN
00483                   f=fade(angle/Outflow%open_angle, .1d0)
00484                   vddr=f*velocity/Outflow%thickness
00485                   dq(ivx:ivy)=Outflow%density*(vddr*(r-Outflow%r_inner)*rotate_z2D((/sign(cos(angle),rpos(1)), sin(angle)/),Outflow%phi+precessPhase)+Outflow%source_vel(1:2))-q(ivx:ivy) !rotate coordinates again so that vector along Outflow-axis is along x-axis
00486                ELSE
00487                   dq(ivx:ivy)=Outflow%density*(Outflow%source_vel(1:2))-q(ivx:ivy)
00488                END IF
00489             END IF
00490          END IF
00491 
00492       CASE(3)
00493 
00494          rpos=rotate_z(rotate_x(pos,-Outflow%phi),-Outflow%theta) !rotate coordinates again so that vector along Outflow-axis is along x-axis
00495          rpos=rotate_z(rotate_x(rpos,-PrecessPhase),-Outflow%PrecessAmplitude)
00496          x=Outflow%offset+abs(rpos(1)) !distance to 'center' for Outflows
00497          s2=rpos(2)**2+rpos(3)**2
00498          s=sqrt(s2)
00499 
00500          IF (Outflow%open_angle == 0) THEN
00501             IF (s < Outflow%radius .AND. x < Outflow%r_outer) THEN
00502                IF (Outflow%B > 0) THEN
00503                   IF (.NOT. ALL(Outflow%source_vel==0)) THEN
00504                      write(*,*) "moving magnetized Outflows not supported"
00505                      STOP
00506                   END IF
00507                   !First calculate existing b field in terms of local rotated cylindrical coordinates
00508                   !                b_curr=rotate_z(rotate_x(q(iBx:iBz),-Outflow%phi),-Outflow%theta)
00509                   !                b_curr(1)=0d0 !We don't want to change B along the axis
00510                   IF (s <= .8*Outflow%radius) THEN !magnetic radius inside Outflow beam
00511                      mag_press=half*(Outflow%B*s/(.8*Outflow%radius))**2
00512                      !                   dq(iBx:iBz)=rotate_x(rotate_z(Outflow%B/(.8*Outflow%radius)*(/0d0,-pos(3),pos(2)/)-b_curr,Outflow%theta),Outflow%phi)
00513                   ELSE
00514                      mag_press=half*(Outflow%B*(Outflow%radius-s)/(.2*Outflow%radius))**2
00515                      !                   dq(iBx:iBz)=rotate_x(rotate_z(Outflow%B*(Outflow%radius-s)/(.2*Outflow%radius*s)*(/0d0,-pos(3),pos(2)/)-b_curr,Outflow%theta),Outflow%phi)
00516                   END IF
00517 
00518                ELSE
00519                   mag_press=0
00520                END IF
00521                pT=Outflow%density*Outflow%temperature
00522                therm_press=pT-mag_press
00523                rho_Outflow=therm_press/Outflow%temperature  
00524                dq(1)=rho_Outflow-q(1)
00525                IF (iE .ne. 0) dq(iE)=therm_press-q(iE)
00526                f=fade(s/Outflow%radius, .1d0)
00527                vddr=f*velocity/Outflow%thickness
00528                IF (x > Outflow%buffer) THEN
00529                   dq(ivx:ivz)=rho_Outflow*(vddr*(x-Outflow%buffer)*rotate_x(rotate_z(rotate_x(rotate_z((/sign(1d0,rpos(1)), 0d0,0d0/),Outflow%theta),Outflow%phi),Outflow%PrecessAmplitude),PrecessPhase)+Outflow%source_vel)-q(ivx:ivz) !rotate coordinates again so that vector along Outflow-axis is along x-axis
00530                ELSE
00531                   dq(ivx:ivz)=rho_Outflow*Outflow%source_vel-q(ivx:ivz)
00532                END IF
00533             END IF
00534          ELSE
00535             angle=atan(s/x)
00536             r=sqrt(s2+x**2)
00537             IF (angle < Outflow%open_angle .AND. r < Outflow%r_outer) THEN
00538                dq(1)=Outflow%density-q(1)
00539                IF (iE .ne. 0) dq(iE)=Outflow%density*Outflow%temperature-q(iE)
00540                IF (r > Outflow%r_inner) THEN
00541                   f=fade(angle/Outflow%open_angle, .1d0)
00542                   vddr=f*velocity/Outflow%thickness
00543                   dq(ivx:ivz)=Outflow%density*(vddr*(r-Outflow%r_inner)*rotate_x(rotate_z(rotate_x(rotate_z((/sign(cos(angle),rpos(1)), sin(angle)*rpos(2)/s, sin(angle)*rpos(3)/s/),Outflow%theta),Outflow%phi),Outflow%PrecessAmplitude),PrecessPhase)+Outflow%source_vel)-q(ivx:ivz) !rotate coordinates again so that vector along Outflow-axis is along x-axis
00544                ELSE
00545                   dq(ivx:ivz)=Outflow%density*Outflow%source_vel-q(ivx:ivz)
00546                END IF
00547             END IF
00548          END IF
00549       END SELECT
00550    END SUBROUTINE OutflowLaunch
00551 
00552    SUBROUTINE OutflowSource(Outflow, pos, q, dq, t, dt)
00553       REAL(KIND=qPrec), DIMENSION(:) :: pos,q, dq
00554       REAL(KIND=qPrec) :: s, r, angle, w(3), vddr, pr, pT, mag_press,therm_press,rz,rz2,s2,x,t,dt,rho_Outflow,f,b_curr(3),precessPhase,oscPhase,velocity,rpos(3)
00555       TYPE(OutflowDef) :: Outflow
00556       pos=pos-Outflow%position !refere to the relative position of source
00557       dq=0
00558       precessphase=0
00559       oscphase=0
00560       IF (Outflow%precessperiod > 0) PrecessPhase=Outflow%PrecessPhase+2*Pi*(t-Outflow%begin)/Outflow%precessPeriod
00561       IF (Outflow%oscPeriod > 0) oscPhase=Outflow%oscPhase+2d0*Pi*(t-Outflow%begin)/Outflow%oscPeriod
00562       velocity=Outflow%velocity+sin(oscPhase)*Outflow%oscAmplitude
00563       !    write(*,*) t,dt,precessphase,oscphase,velocity
00564       SELECT CASE(nDim)
00565       CASE(2)
00566 
00567          rpos(1:2)=rotate_z2D(pos(1:2),-Outflow%phi-PrecessPhase) !rotate to coordinates in Outflow frame
00568          x=Outflow%offset+abs(rpos(1)) !distance to 'center' for Outflows
00569          s2=rpos(2)**2
00570          s=abs(rpos(2))
00571 
00572          IF (Outflow%open_angle == 0) THEN
00573             IF (s < Outflow%radius .AND. x < Outflow%r_outer) THEN
00574                IF (Outflow%B > 0) THEN
00575                   write(*,*) "Magnetized Outflows in 2D not supported"
00576                   STOP
00577                END IF
00578                dq(1)=Outflow%density-q(1)
00579                therm_press=Outflow%density*Outflow%temperature
00580                IF (iE .ne. 0) dq(iE)=therm_press-q(iE)
00581                f=fade(s/Outflow%radius, .5d0)
00582                vddr=f*velocity/Outflow%thickness
00583                IF (x > Outflow%buffer) THEN
00584                   dq(ivx:ivy)=Outflow%density*(vddr*(x-Outflow%buffer)*rotate_z2D((/sign(1d0,rpos(1)), 0d0/),Outflow%phi+precessPhase)+Outflow%source_vel(1:2))-q(ivx:ivy) !rotate coordinates again so that vector along Outflow-axis is along x-axis
00585                ELSE
00586                   dq(ivx:ivy)=Outflow%density*(Outflow%source_vel(1:2))-q(ivx:ivy)
00587                END IF
00588             END IF
00589          ELSEIF (abs(Outflow%open_angle-half*Pi) < .01) THEN
00590             !          angle=atan(s/x)
00591             r=sqrt(s2+x**2)
00592             IF (r < Outflow%r_outer) THEN
00593                dq(1)=Outflow%density-q(1)
00594                IF (iE .ne. 0) dq(iE)=Outflow%density*Outflow%temperature-q(iE)
00595                IF (r > Outflow%r_inner) THEN
00596                   vddr=velocity/Outflow%thickness
00597                   dq(ivx:ivy)=Outflow%density*(vddr*(r-Outflow%r_inner)*rotate_z2D(rpos(1:2)/r,Outflow%phi+precessPhase)+Outflow%source_vel(1:2))-q(ivx:ivy) !rotate coordinates again so that vector along Outflow-axis is along x-axis
00598                ELSE
00599                   dq(ivx:ivy)=Outflow%density*(Outflow%source_vel(1:2))-q(ivx:ivy)
00600                END IF
00601             END IF
00602          ELSE
00603             angle=atan(s/x)
00604             r=sqrt(s2+x**2)
00605             IF (angle < Outflow%open_angle .AND. r < Outflow%r_outer) THEN
00606                dq(1)=Outflow%density-q(1)
00607                IF (iE .ne. 0) dq(iE)=Outflow%density*Outflow%temperature-q(iE)
00608                IF (r > Outflow%r_inner) THEN
00609                   f=fade(angle/Outflow%open_angle, .1d0)
00610                   vddr=f*velocity/Outflow%thickness
00611                   dq(ivx:ivy)=Outflow%density*(vddr*(r-Outflow%r_inner)*rotate_z2D((/sign(cos(angle),rpos(1)), sin(angle)/),Outflow%phi+precessPhase)+Outflow%source_vel(1:2))-q(ivx:ivy) !rotate coordinates again so that vector along Outflow-axis is along x-axis
00612                ELSE
00613                   dq(ivx:ivy)=Outflow%density*(Outflow%source_vel(1:2))-q(ivx:ivy)
00614                END IF
00615             END IF
00616          END IF
00617 
00618       CASE(3)
00619 
00620          rpos=rotate_z(rotate_x(pos,-Outflow%phi),-Outflow%theta) !rotate coordinates again so that vector along Outflow-axis is along x-axis
00621          rpos=rotate_z(rotate_x(rpos,-PrecessPhase),-Outflow%PrecessAmplitude)
00622          x=Outflow%offset+abs(rpos(1)) !distance to 'center' for Outflows
00623          s2=rpos(2)**2+rpos(3)**2
00624          s=sqrt(s2)
00625 
00626          IF (Outflow%open_angle == 0) THEN
00627             IF (s < Outflow%radius .AND. x < Outflow%r_outer) THEN
00628                IF (Outflow%B > 0) THEN
00629                   IF (.NOT. ALL(Outflow%source_vel==0)) THEN
00630                      write(*,*) "moving magnetized Outflows not supported"
00631                      STOP
00632                   END IF
00633                   !First calculate existing b field in terms of local rotated cylindrical coordinates
00634                   !                b_curr=rotate_z(rotate_x(q(iBx:iBz),-Outflow%phi),-Outflow%theta)
00635                   !                b_curr(1)=0d0 !We don't want to change B along the axis
00636                   IF (s <= .8*Outflow%radius) THEN !magnetic radius inside Outflow beam
00637                      mag_press=half*(Outflow%B*s/(.8*Outflow%radius))**2
00638                      !                   dq(iBx:iBz)=rotate_x(rotate_z(Outflow%B/(.8*Outflow%radius)*(/0d0,-pos(3),pos(2)/)-b_curr,Outflow%theta),Outflow%phi)
00639                   ELSE
00640                      mag_press=half*(Outflow%B*(Outflow%radius-s)/(.2*Outflow%radius))**2
00641                      !                   dq(iBx:iBz)=rotate_x(rotate_z(Outflow%B*(Outflow%radius-s)/(.2*Outflow%radius*s)*(/0d0,-pos(3),pos(2)/)-b_curr,Outflow%theta),Outflow%phi)
00642                   END IF
00643 
00644                ELSE
00645                   mag_press=0
00646                END IF
00647                pT=Outflow%density*Outflow%temperature
00648                therm_press=pT-mag_press
00649                rho_Outflow=therm_press/Outflow%temperature  
00650                dq(1)=rho_Outflow-q(1)
00651                IF (iE .ne. 0) dq(iE)=therm_press-q(iE)
00652                f=fade(s/Outflow%radius, .1d0)
00653                vddr=f*velocity/Outflow%thickness
00654                IF (x > Outflow%buffer) THEN
00655                   dq(ivx:ivz)=rho_Outflow*(vddr*(x-Outflow%buffer)*rotate_x(rotate_z(rotate_x(rotate_z((/sign(1d0,rpos(1)), 0d0,0d0/),Outflow%theta),Outflow%phi),Outflow%PrecessAmplitude),PrecessPhase)+Outflow%source_vel)-q(ivx:ivz) !rotate coordinates again so that vector along Outflow-axis is along x-axis
00656                ELSE
00657                   dq(ivx:ivz)=rho_Outflow*Outflow%source_vel-q(ivx:ivz)
00658                END IF
00659             END IF
00660          ELSE
00661             angle=atan(s/x)
00662             r=sqrt(s2+x**2)
00663             IF (angle < Outflow%open_angle .AND. r < Outflow%r_outer) THEN
00664                dq(1)=Outflow%density-q(1)
00665                IF (iE .ne. 0) dq(iE)=Outflow%density*Outflow%temperature-q(iE)
00666                IF (r > Outflow%r_inner) THEN
00667                   f=fade(angle/Outflow%open_angle, .1d0)
00668                   vddr=f*velocity/Outflow%thickness
00669                   dq(ivx:ivz)=Outflow%density*(vddr*(r-Outflow%r_inner)*rotate_x(rotate_z(rotate_x(rotate_z((/sign(cos(angle),rpos(1)), sin(angle)*rpos(2)/s, sin(angle)*rpos(3)/s/),Outflow%theta),Outflow%phi),Outflow%PrecessAmplitude),PrecessPhase)+Outflow%source_vel)-q(ivx:ivz) !rotate coordinates again so that vector along Outflow-axis is along x-axis
00670                ELSE
00671                   dq(ivx:ivz)=Outflow%density*Outflow%source_vel-q(ivx:ivz)
00672                END IF
00673             END IF
00674          END IF
00675       END SELECT
00676    END SUBROUTINE OutflowSource
00677 
00678    FUNCTION emf_Launch_2D(Outflow,pos)
00679       REAL(KIND=qPREC) :: pos(:),tor,s2,s,x,t,precessPhase,oscPhase
00680       TYPE(OutflowDef) :: Outflow
00681       REAL(KIND=qPREC) :: emf_Launch_2D
00682       emf_launch_2d=0
00683       print*, '2D emfs not implemented for outflows'
00684       STOP
00685    END FUNCTION emf_Launch_2D
00686 
00687    FUNCTION emf_Launch_3D(Outflow,pos)
00688       REAL(KIND=qPREC) :: pos(:),tor,s2,s,x,t,precessPhase,oscPhase
00689       TYPE(OutflowDef) :: Outflow
00690       REAL(KIND=qPREC) :: emf_Launch_3D(3)
00691       pos=pos-Outflow%position !refere to the relative position of source
00692       precessphase=0
00693       oscphase=0
00694       IF (Outflow%precessperiod > 0) PrecessPhase=Outflow%PrecessPhase
00695       IF (Outflow%oscPeriod > 0) oscPhase=Outflow%oscPhase
00696       pos=rotate_z(rotate_x(pos,-Outflow%phi),-Outflow%theta) !rotate coordinates again so that vector along Outflow-axis is along x-axis
00697       pos=rotate_z(rotate_x(pos,-precessPhase),-Outflow%PrecessAmplitude)
00698       !    velocity=Outflow%velocity+sin(oscPhase)*Outflow%oscAmplitude
00699 
00700       x=Outflow%offset+abs(pos(1)) !distance to 'center' for Outflows
00701       s2=pos(2)**2+pos(3)**2
00702       s=sqrt(s2)
00703       emf_launch_3D=0
00704       IF (x < Outflow%r_outer) THEN
00705          IF (s <= .8*Outflow%radius) THEN
00706             tor=Outflow%B*(.8*Outflow%radius**2-s2)/(2d0*.8*Outflow%radius)
00707             emf_launch_3D=tor*rotate_x(rotate_z(rotate_x(rotate_z((/1d0,0d0,0d0/),Outflow%theta),Outflow%phi),Outflow%PrecessAmplitude),precessPhase)
00708          ELSEIF (s < Outflow%radius) THEN
00709             tor=Outflow%B*(Outflow%radius-s)**2/(2d0*(.2)*Outflow%radius)
00710             emf_launch_3D=tor*rotate_x(rotate_z(rotate_x(rotate_z((/1d0,0d0,0d0/),Outflow%theta),Outflow%phi),Outflow%PrecessAmplitude),PrecessPhase)
00711          END IF
00712       END IF
00713    END FUNCTION emf_Launch_3D
00714 
00715 
00716    
00717   FUNCTION emf_source_3D(Outflow,pos,t,dt)
00718       REAL(KIND=qPREC) :: pos(:),tor,s2,s,x,t,dt,vddr,PrecessPhase,OscPhase,velocity
00719       TYPE(OutflowDef) :: Outflow
00720       REAL(KIND=qPREC), DIMENSION(3) :: emf_source_3D
00721       pos=pos-Outflow%position !refere to the relative position of source
00722       precessphase=0
00723       oscphase=0
00724       IF (Outflow%precessperiod > 0) PrecessPhase=Outflow%PrecessPhase+2*Pi*(t-Outflow%begin)/Outflow%precessPeriod
00725       IF (Outflow%oscPeriod > 0) oscPhase=Outflow%oscPhase+2d0*Pi*(t-Outflow%begin)/Outflow%oscPeriod
00726       pos=rotate_z(rotate_x(pos,-Outflow%phi),-Outflow%theta) !rotate coordinates again so that vector along Outflow-axis is along x-axis
00727       pos=rotate_z(rotate_x(pos,-precessPhase),-Outflow%PrecessAmplitude)
00728       velocity=Outflow%velocity+sin(oscPhase)*Outflow%oscAmplitude
00729 
00730       x=Outflow%offset+abs(pos(1)) !distance to 'center' for Outflows
00731       s2=pos(2)**2+pos(3)**2
00732       s=sqrt(s2)
00733       emf_source_3D=0
00734       IF (x < Outflow%r_outer .AND. x > Outflow%buffer) THEN
00735          vddr=velocity/Outflow%thickness
00736          IF (s <= .8*Outflow%radius) THEN
00737             tor=vddr*Outflow%B*(.8*Outflow%radius**2-s2)/(2d0*.8*Outflow%radius)
00738             emf_source_3D=tor*rotate_x(rotate_z(rotate_x(rotate_z((/1d0,0d0,0d0/),Outflow%theta),Outflow%phi),Outflow%PrecessAmplitude),PrecessPhase)*dt
00739          ELSEIF (s < Outflow%radius) THEN
00740             tor=vddr*Outflow%B*(Outflow%radius-s)**2/(2d0*(.2)*Outflow%radius)
00741             emf_source_3D=tor*rotate_x(rotate_z(rotate_x(rotate_z((/1d0,0d0,0d0/),Outflow%theta),Outflow%phi),Outflow%PrecessAmplitude),PrecessPhase)*dt
00742          END IF
00743       END IF
00744    END FUNCTION emf_source_3D
00745 
00746   FUNCTION emf_source_2D(Outflow,pos,t,dt)
00747       REAL(KIND=qPREC) :: pos(:),tor,s2,s,x,t,dt,vddr,PrecessPhase,OscPhase,velocity
00748       TYPE(OutflowDef) :: Outflow
00749       REAL(KIND=qPREC) :: emf_source_2D
00750       pos=pos-Outflow%position !refere to the relative position of source
00751       emf_source_2d=0
00752       print*, 'emf source not supported for outflows in 2d'
00753       STOP
00754    END FUNCTION emf_source_2D
00755 
00756 
00759    SUBROUTINE SetOutflowGeometry(Outflow)
00760       TYPE(OutflowDef) :: Outflow
00761       !      Outflow%temperature=Outflow%temperature/EosConstants
00762       IF (Outflow%open_angle == 0) THEN
00763          Outflow%offset=0
00764          Outflow%r_inner=Outflow%radius
00765          Outflow%r_outer=Outflow%buffer+Outflow%thickness
00766       ELSEIF (abs(Outflow%open_angle-half*Pi) < .01) THEN
00767          Outflow%offset=0
00768          Outflow%buffer=0
00769          Outflow%r_inner=Outflow%radius
00770          Outflow%r_outer=Outflow%radius+Outflow%thickness
00771       ELSE
00772          Outflow%offset=max(0d0,Outflow%radius/tan(Outflow%open_angle)-Outflow%buffer)
00773          Outflow%buffer=min(Outflow%buffer,Outflow%offset)
00774          Outflow%r_inner=sqrt((Outflow%offset+Outflow%buffer)**2+Outflow%radius**2)
00775          Outflow%r_outer=Outflow%r_inner+Outflow%thickness
00776       END IF
00777       IF (Outflow%open_angle==0) THEN
00778          Outflow%radius_needed=sqrt(Outflow%r_inner**2+Outflow%r_outer**2)
00779       ELSE
00780          Outflow%radius_needed=sqrt((Outflow%r_outer*sin(Outflow%open_angle))**2+(Outflow%r_outer*cos(Outflow%open_angle)-Outflow%offset)**2)
00781       END IF
00782    END SUBROUTINE SetOutflowGeometry
00783 
00784    SUBROUTINE SetOutflowBounds(Outflow)
00785       TYPE(OutflowDef) :: Outflow
00786       Outflow%xBounds=0
00787       Outflow%xBounds(1:nDim,1)=Outflow%position(1:nDim)-Outflow%radius_needed
00788       Outflow%xBounds(1:nDim,2)=Outflow%position(1:nDim)+Outflow%radius_needed
00789    END SUBROUTINE SetOutflowBounds
00790 
00791 END MODULE Outflows
00792 
00793 
 All Classes Files Functions Variables