musrsim/accel/src/SUB_INTEGR_6.FOR

1008 lines
29 KiB
Fortran

OPTIONS /EXTEND_SOURCE
SUBROUTINE READ_INFOS
c =====================
logical map_error /.false./
COMMON /map_error/ map_error
CALL READ_INFO_6 ! Die Reihenfolge der 'READ_INFO_x'-calls
CALL READ_INFO_5 ! ist wichtig (wegen Uebermittlung der
CALL READ_INFO_4 ! Uebergabebereiche)
CALL READ_INFO_3
CALL READ_INFO_2
CALL READ_INFO_1
if (map_error) then
write(*,*) '-----------------------------------------------------------'
write(*,*)
write(*,*) ' => Files ''MAP_DEF_x.INC'' (x == Nr. der Mappe) im Directory'
write(*,*) ' ''accel$sourcedirectory:'' angleichen und anschliessend'
write(*,*) ' das Programm mit ''LINKACV'' am DCL Prompt neu kompilieren'
write(*,*) ' und linken.'
write(*,*)
write(*,*) '-----------------------------------------------------------'
write(*,*)
STOP
endif
END
c===============================================================================
OPTIONS /EXTEND_SOURCE
SUBROUTINE READ_INFO_6
c ======================
IMPLICIT NONE
character*1 Nr
parameter (Nr='6')
INCLUDE 'accel$sourcedirectory:COM_GEO.INC'
INCLUDE 'accel$sourcedirectory:COM_HVs.INC'
INCLUDE 'accel$sourcedirectory:COM_ACCEL.INC'
INCLUDE 'accel$sourcedirectory:COM_DIRS.INC'
INCLUDE 'accel$sourcedirectory:MAP_DEF_6.INC'
c - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
c nur um Mappeninfos fuer Runs vor Run 10 weiterhin lesen zu koennen
c (wurden ab Version 1.2.1 ersetzt durch outerDy_TgtHolder, outerDz_TgtHolder):
real Dy_TgtHolder, Dz_TgtHolder
c - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
namelist /voltages/ UTgt, UGua, UGrid1, UG1, UG2, zero
namelist /geometry/
+ Dy_Foil,Dz_Foil,xEnd_TgtHolder,Dy_TgtHolder,Dz_TgtHolder,
+ outerDy_TgtHolder,outerDz_TgtHolder,
+ innerDy1_TgtHolder,innerDz1_TgtHolder,innerDy2_TgtHolder,innerDz2_TgtHolder,
+ xStart_Guardring,xEnd_Guardring,innerDy_Guardring,outerDy_Guardring,
+ innerDz_Guardring,outerDz_Guardring,
+ xPosition_Grid1,distance_wires1,dWires1,y_Pos_firstWire1,y_Pos_lastWire1,
+ xStart_Balken,xEnd_Balken,Dy_Balken,
+ innerDz_Balken,outerDz_Balken,
+ xStart_Gridframe1,xEnd_Gridframe1,innerDy_Gridframe1,outerDy_Gridframe1,
+ innerDz_Gridframe1,outerDz_Gridframe1,
+ xPosition_Grid2,distance_wires2,dWires2,y_Pos_firstWire2,y_Pos_lastWire2,
+ xStart_Gridframe2,xEnd_Gridframe2,innerDy_Gridframe2,outerDy_Gridframe2,
+ innerDz_Gridframe2,outerDz_Gridframe2,
+ xHeShield,rHeShield
c - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
logical map_error
COMMON /map_error/ map_error
c the grid characteristics (as read from the INFO-file):
real Dx,Dy,Dz
real x_iEQ1, ymax,zmax ! xmax wird in MAP_DEF_6.INC definiert
namelist /grid_info/
+ Dx,Dy,Dz, imax,jmax,kmax, x_iEQ1, xmin,xmax,ymax,zmax
integer iostat
c - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
c Versuchen, ein Mappenfile mit fixer HV-Einstellung zu oeffnen (d.h. ein
c Mappenfile ohne Erweiterung '_Tgt' bzw. '_Gua' bzw. '_Gi1'). Gelingt dies, so
c verwende fuer die Hochspannungen an Target und 1. Gitter die Spannungen,
c fuer die die Mappe gerechnet wurde (-> 'HVs_from_map = .true.). Ansonsten
c verwende die von ACCEL.INPUT eingelesenen Vorgaben fuer die HV-Schleifen:
open (lunREAD,file=mappenName//'_'//Nr,defaultfile=mappenDir//':.MAPPE',
+ readonly,status='old',iostat=iostat)
if (iostat.EQ.0) then
HVs_from_map = .true. ! Default = .false.
close (lunREAD)
endif
c Einlesen der Mappen-Informationen (bei Mappe 6 auch Spannungen und Geometrie
c des Targetbereiches):
open (lunREAD,file=mappenName//'_'//Nr,defaultfile=mappenDir//':.INFO',
+ readonly,status='old')
read(lunREAD,nml=grid_info)
rewind (lunREAD)
UGua = -1.E10
read(lunREAD,nml=voltages)
if (UGua.NE.-1.E10) then
! Es wurde ein Spannungswert fuer den Guardring angegeben.
! => Guardring hat eigene Spannungseinstellung
freeGuard = .true.
else
freeGuard = .false.
par(1,UGuard) = 0
par(2,UGuard) = 0
par(3,UGuard) = 1
endif
rewind (lunREAD)
read(lunREAD,nml=geometry)
close (lunREAD)
if (outerDy_TgtHolder.EQ.-1.E10) outerDy_TgtHolder = Dy_TgtHolder
if (outerDz_TgtHolder.EQ.-1.E10) outerDz_TgtHolder = Dz_TgtHolder
if (innerDy1_TgtHolder.EQ.-1.E10) innerDy1_TgtHolder = Dy_Foil
if (innerDz1_TgtHolder.EQ.-1.E10) innerDz1_TgtHolder = Dz_Foil
if (innerDy2_TgtHolder.EQ.-1.E10) innerDy2_TgtHolder = Dy_Foil
if (innerDz2_TgtHolder.EQ.-1.E10) innerDz2_TgtHolder = Dz_Foil
if (HVs_from_map) then
! Die eingelesenen Spannungen an die Variablen fuer die
! zugehoerigen Schleifen weitergeben (nur bei Mappe 6):
par(1,UTarget) = UTgt
par(2,UTarget) = UTgt
par(3,UTarget) = 1
if (freeGuard) then
par(1,UGuard) = UGua
par(2,UGuard) = UGua
par(3,UGuard) = 1
endif
par(1,UGi1) = UG1
par(2,UGi1) = UG1
par(3,UGi1) = 1
endif
c eingelesene imax, jmax und kmax um 1 reduzieren, da in 'ACCEL' die Feldindizes
c ab 0 laufen, bei 'RELAX3D' jedoch ab 1:
imax = imax-1
jmax = jmax-1
kmax = kmax-1
c Umrechnen der Koordinaten, wie sie von 'BESCHL-INIT' ('RELAX3D') verwendet
c werden (Ursprung in Targetfolienmitte) in System mit Ursprung auf der Kryo-Achse:
xFoil = -xHeShield
xHeShield = 0.
xmin = xmin + xFoil
xmax = xmax + xFoil
xStartMap6 = xmin
xEndMap6 = xmax
xEnd_TgtHolder = xEnd_TgtHolder + xFoil
xStart_Guardring = xStart_Guardring + xFoil
xEnd_Guardring = xEnd_Guardring + xFoil
xPosition_Grid1 = xPosition_Grid1 + xFoil
xStart_Gridframe1 = xStart_Gridframe1 + xFoil
xEnd_Gridframe1 = xEnd_Gridframe1 + xFoil
xStart_Balken = xStart_Balken + xFoil
xEnd_Balken = xEnd_Balken + xFoil
xPosition_Grid2 = xPosition_Grid2 + xFoil
xStart_Gridframe2 = xStart_Gridframe2 + xFoil
xEnd_Gridframe2 = xEnd_Gridframe2 + xFoil
rWires1 = dWires1/2.
rQuadWires1 = rWires1*rWires1
rWires2 = dWires2/2.
rQuadWires2 = rWires2*rWires2
C DER FOLGENDE ABSCHNITT WURDE HERAUSKOMMENTIERT, DA ES MITTLERWEILE VERSCHIEDEN
C GROSSE POTENTIALMAPPEN GIBT UND DIE MAPPENDIMENSIONEN DAHER SOWIESO VARIABEL
C GEHALTEN WERDEN MUESSEN. DIE VERWENDUNG VON PARAMETERN IST LEIDER NICHT
C MEHR MOEGLICH. ('LEIDER' WEGEN DER ERHOEHTEN RECHENZEIT):
C
Cc checken, ob die Charakteristika der einzulesenden Mappe mit den Vorgaben der
Cc Integrationsroutinen uebereinstimmen:
Cc (Die Laenge von Mappe 6 ist abhaengig von der Position des Moderators relativ
Cc zur Kryoachse.)
C
C if (
Cc + imax.NE.imax_ .OR.
C + jmax.NE.jmax_ .OR. kmax.NE.kmax_ .OR.
C + Dx.NE.Dx_ .OR. Dy.NE.Dy_ .OR. Dz.NE.Dz_
Cc + .OR. xmin.NE.xmin_
C + ) then
C write(*,*) '-----------------------------------------------------------'
C if (.NOT.map_error) then
C write(*,*) ' Feldgroessen der eingelesenen Mappe und des reservierten'
C write(*,*) ' Speichers stimmen nicht ueberein:'
C write(*,*)
C endif
C write(*,*) ' MAPPE '//Nr//': '//mappenName(1:nameLength)//'_'//Nr//'.MAPPE'
C write(*,*) ' Mappe: imax, jmax ,kmax = ',imax ,jmax ,kmax
C write(*,*) ' Dx ,Dy ,Dz = ',Dx ,Dy ,Dz
C write(*,*) ' Speicher: imax_,jmax_,kmax_ = ',imax ,jmax_,kmax_
C write(*,*) ' Dx_ ,Dy_ ,Dz_ = ',Dx_ ,Dy_ ,Dz_
C write(*,*)
C map_error = .true.
C endif
C
C if (map_error) RETURN ! kann auch in anderem 'READ_MAP_x' gesetzt worden sein
c checken, ob der reservierte Speicherplatz ausreicht:
if ((imax+1)*(jmax+1)*(kmax+1).GT.maxmaxmem+1) then
write(*,*)
write(*,*) 'reservierter Speicher ist nicht ausreichend fuer Mappe',Nr
write(*,*)
write(*,*) '(imax+1)*(jmax+1)*(kmax+1) = ',(imax+1)*(jmax+1)*(kmax+1)
write(*,*) 'maxmaxmem + 1 = ',maxmaxmem + 1
write(*,*)
write(*,*) '=> ''maxmaxmem'' in accel$sourcedirectory:MAPMAP.INC angleichen'
write(*,*) ' und Programm mit ''LINKACV'' am DCL-Prompt neu kompilieren'
write(*,*) ' und linken.'
write(*,*)
call exit
endif
END
c===============================================================================
OPTIONS /EXTEND_SOURCE
SUBROUTINE READ_MAP_6
c =====================
IMPLICIT NONE
character*1 Nr
parameter (Nr='6')
INCLUDE 'accel$sourcedirectory:MAP_DEF_6.INC'
INCLUDE 'accel$sourcedirectory:COM_DIRS.INC'
INCLUDE 'accel$sourcedirectory:COM_ACCEL.INC'
integer i,j,k, ihelp, iostat
c Einlesen der Potentialmappe:
open (lunRead,file=mappenName//'_'//Nr,
+ defaultfile=mappenDir//':.MAPPE',status='old',
+ form='unformatted',recl=imax+1,readonly)
c + form='unformatted',recl=imax,readonly)
write(*,*) 'reading '//mappenName(1:nameLength)//'_'//Nr//'.MAPPE ...'
do k = 0, kmax
do j = 0, jmax
c read(lunREAD,iostat=iostat) (map(i,j,k),i=0,imax)
ihelp = (k*(jmax+1)+j)*(imax+1)
read(lunREAD,iostat=iostat) (map(ihelp+i),i=0,imax)
if (iostat.NE.0) then
write(*,*)
write(*,999) i,j,k,iostat
STOP
endif
enddo
enddo
close(lunREAD)
999 format(x/'error reading grid point (i,j,k) = ('i4','i4','
+ i4')'/'iostat = 'i4/)
c da die Anodenbereiche bei RELAX3D negativ kodiert sind, nimm die
c Absolutbetraege:
ihelp = 0
do k=0, kmax
do j=0, jmax
do i=0, imax
c map(i,j,k) = abs(map(i,j,k))
map(ihelp) = abs(map(ihelp))
ihelp = ihelp + 1
enddo
enddo
enddo
RETURN
END
c===============================================================================
OPTIONS /EXTEND_SOURCE
SUBROUTINE ADD_MAP_6
c ====================
IMPLICIT NONE
character*1 Nr
parameter (Nr='6')
INCLUDE 'accel$sourcedirectory:MAP_DEF_6.INC'
INCLUDE 'accel$sourcedirectory:COM_HVs.INC'
INCLUDE 'accel$sourcedirectory:COM_DIRS.INC'
INCLUDE 'accel$sourcedirectory:COM_ACCEL.INC'
real read_memory(0:100)
COMMON /read_memory/ read_memory ! COMMON nur, damit nicht jede
! Mappe extra Speicher belegt.
integer i,j,k, ihelp, iostat
c Einlesen der '_Tgt_Nr'-Potentialmappe:
if (mappenName.EQ.'RUN9') then
open (lunRead,file='RUN6_NEW_Tgt_'//Nr,
+ defaultfile=mappenDir//':.MAPPE',status='old',
+ form='unformatted',recl=imax+1,readonly)
c + form='unformatted',recl=imax,readonly)
else
open (lunRead,file=mappenName//'_Tgt_'//Nr,
+ defaultfile=mappenDir//':.MAPPE',status='old',
+ form='unformatted',recl=imax+1,readonly)
c + form='unformatted',recl=imax,readonly)
endif
write(*,*) 'constructing map '//Nr//' ...'
do k = 0, kmax
do j = 0, jmax
c read(lunREAD,iostat=iostat) (map(i,j,k),i=0,imax)
ihelp = (k*(jmax+1)+j)*(imax+1)
read(lunREAD,iostat=iostat) (map(ihelp+i),i=0,imax)
if (iostat.NE.0) then
write(*,*)
write(*,999) i,j,k,iostat
STOP
endif
enddo
enddo
close(lunREAD)
999 format(x/'error reading grid point (i,j,k) = ('i4','i4','
+ i4')'/'iostat = 'i4/)
c Angleichen der Potentialmappe an UTgt:
ihelp = 0
do k=0, kmax
do j=0, jmax
do i=0, imax
c map(i,j,k) = UTgt*abs(map(i,j,k))
map(ihelp) = UTgt*abs(map(ihelp))
ihelp = ihelp + 1
enddo
enddo
enddo
c Einlesen der '_Gua_Nr'-Potentialmappe und Angleichen an UGua:
if (freeGuard) then
open (lunRead,file=mappenName//'_Gua_'//Nr,
+ defaultfile=mappenDir//':.MAPPE',status='old',
+ form='unformatted',recl=imax+1,readonly)
c + form='unformatted',recl=imax,readonly)
ihelp = 0
do k = 0, kmax
do j = 0, jmax
read(lunRead,iostat=iostat) (read_memory(i),i=0,imax)
if (iostat.NE.0) then
write(*,*)
write(*,999) i,j,k,iostat
STOP
endif
do i=0, imax
map(ihelp) = map(ihelp) + UGua*abs(read_memory(i))
ihelp = ihelp + 1
enddo
enddo
enddo
close(lunREAD)
endif
c Einlesen der '_Gi1_Nr'-Potentialmappe und Angleichen an UG1:
open (lunRead,file=mappenName//'_Gi1_'//Nr,
+ defaultfile=mappenDir//':.MAPPE',status='old',
+ form='unformatted',recl=imax+1,readonly)
c + form='unformatted',recl=imax,readonly)
ihelp = 0
do k = 0, kmax
do j = 0, jmax
read(lunRead,iostat=iostat) (read_memory(i),i=0,imax)
if (iostat.NE.0) then
write(*,*)
write(*,999) i,j,k,iostat
STOP
endif
do i=0, imax
map(ihelp) = map(ihelp) + UG1*abs(read_memory(i))
ihelp = ihelp + 1
enddo
enddo
enddo
close(lunREAD)
RETURN
END
c===============================================================================
OPTIONS /EXTEND_SOURCE
SUBROUTINE INTEGRATIONSSTEP_RUNGE_KUTTA_6(dt)
c =============================================
IMPLICIT NONE
SAVE
character*1 Nr
parameter (Nr='6')
c Diese Subroutine berechnet zu einem vorgegebenen Zeitschritt dt den
c Integrationsschritt zweimal: einmal direkt mit dt und einmal ueber zwei
c aufeinanderfolgende Schritte mit dt/2. (die beiden dt/2-Schritte werden
c zuerst ausgefuehrt).
c
c Aus der Differenz der beiden Resultate wird eine Abschaetzung fuer den Fehler
c des dt-Schrittes gewonnen, die dazu verwendet wird zu entscheiden, ob der
c Integrationsschritt mit einem verkuerzten Zeitintervall wiederholt werden
c muss, oder ob das Zeitintervall fuer den folgenden ausgedehnt werden kann.
c
c Die beiden Einzelergebnisse aus dem dt- und den beiden dt/2-Schritten, die
c jeweils ueber Runge-Kutta-Rechnung vierter Ordnung erhalten werden, werden
c zum Schluss noch zusammengenommen, um ein Resultat mit Genauigkeit fuenfter
c Ordnung in dt zu erhalten.
c
c Der ganze Ablauf entspricht den Ausfuehrungen in Kapitel 15.2 der NUMERICAL
c RECIPIES: 'Adaptive Stepsize Control for Runge-Kutta' (vgl. Referenz im
c fileheader von 'ACCEL.FOR')
INCLUDE 'accel$sourcedirectory:COM_ACCEL.INC'
INCLUDE 'accel$sourcedirectory:MAP_DEF_6.INC'
integer i ! Zaehlvariable
real dt,dt_half ! zeitl. Aenderung, halbe zeitl. Aenderung
real EFeld0(3), EFeld1(3) ! elektr. Felder
real x1(3),Dx1(3),Dx2(3) ! fuer Ortsintegration
real v1(3),Dv1(3),Dv2(3) ! fuer Geschw.Integration
real xDifferenz(3), vDifferenz(3)
real maxErr_x,maxErr_v,maxErr ! fuer Fehlerbetrachtung
real errCon, safety ! fuer Schrittweitenkontrolle
real pShrink, pGrow ! fuer Schrittweitenkontrolle
PARAMETER (errCon = 6.e-4, safety = .9) ! vgl. Referenz
PARAMETER (pShrink = -.25, pGrow = -.2)
! errCon = (4./safety)**(1/pGrow)
logical flag_dtSmall ! wenn dt kleiner als dtsmall ist und
! der Fehler immer noch zu gross ist.
integer returnCode_EFeld
COMMON /returnCode_EFeld/ returnCode_EFeld
! 1: Testort hinter der Mappe
! 2: Testort neben der Mappe
! 3: Testort vor der Mappe
logical flag_hang
COMMON /flag_hang/ flag_hang
c - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
flag_dtSmall = .false. ! flag resetten
if (dt.lt.dtsmall) dt = dtsmall
c berechne EFeld am aktuellen Ort. Speichere in EFeld0, damit sie wiederverwendet
c werden kann, falls mit kuerzerem Zeitschritt wiederholt werden muss:
call EFeld_6(x,EFeld0,*999)
c...............................................................................
10 continue ! hier gehts wieder von vorne los, falls Zeitschritt dt
! abgeaendert werden muss.
dt_half = dt / 2.
c mache ersten dt/2 - Schritt:
call SINGLESTEP_RUNGE_KUTTA_6(dt_half,EFeld0,x,v, Dx1,Dv1 ,*999)
c berechne EFeld bei x1:
x1(1) = x(1) + Dx1(1)
x1(2) = x(2) + Dx1(2)
x1(3) = x(3) + Dx1(3)
v1(1) = v(1) + Dv1(1)
v1(2) = v(2) + Dv1(2)
v1(3) = v(3) + Dv1(3)
call EFeld_6(x1,EFeld1,*999)
c mache zweiten dt/2 - Schritt:
call SINGLESTEP_RUNGE_KUTTA_6(dt_half,EFeld1,x1,v1, Dx2,Dv2 ,*999)
c Summiere Ergebnisse der beiden dt/2 -Schritte und speichere in Dx1, Dv1:
Dx1(1) = Dx1(1) + Dx2(1)
Dx1(2) = Dx1(2) + Dx2(2)
Dx1(3) = Dx1(3) + Dx2(3)
Dv1(1) = Dv1(1) + Dv2(1)
Dv1(2) = Dv1(2) + Dv2(2)
Dv1(3) = Dv1(3) + Dv2(3)
c mache dt - Schritt:
call SINGLESTEP_RUNGE_KUTTA_6(dt,EFeld0,x,v, Dx2,Dv2 ,*999)
c . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
c Fehlerbetrachtung und gegebenenfalls Berechnung von neuem Ort und neuer
c Geschwindigkeit (falls der Fehler ausserhalb der Toleranz liegt wird Zeit-
c schritt dt verkuerzt und bei Label 10 erneut begonnen):
c Zaehle die Schritte:
steps = Steps + 1
c Fehlerbetrachtung:
c der groesste (absolute bzw. relative) Fehler im Ort soll kleiner als eps_x
c sein, der groesste Fehler in der Geschwindigkeit kleiner als eps_v:
c -> Bestimme den jeweils groessten Fehler der drei Komponenten des Ortes und
c der Geschwindigkeit (dh. die groesste Differenz der Aederungen):
maxErr_x = 0.
maxErr_v = 0.
do i = 1, 3
xDifferenz(i) = Dx1(i)-Dx2(i)
vDifferenz(i) = Dv1(i)-Dv2(i)
if (log_relativ) then
if (Dx1(i).NE.0.) maxErr_x = Max(maxErr_x,Abs(xDifferenz(i)/Dx1(i)))
if (Dv1(i).NE.0.) maxErr_v = Max(maxErr_v,Abs(vDifferenz(i)/Dv1(i)))
else
maxErr_x = Max(maxErr_x,Abs(xDifferenz(i)))
maxErr_v = Max(maxErr_v,Abs(vDifferenz(i)))
endif
enddo
c - Skaliere den jeweils groessten relativen Fehler auf das jeweilige Epsilon:
maxErr_x = maxErr_x / eps_x
maxErr_v = maxErr_v / eps_v
c . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
c der groessere der beiden reskalierten Fehler bestimmt, ob der Integrations-
c schritt mit kleinerem Zeitintervall wiederholt werden muss, bzw. um welchen
c Faktor das Zeitintervall fuer den naechsten Schritt vergroessert werden kann:
c Liegt der Fehler ausserhalb des Toleranzbereiches und ist dt bereits jetzt
c kleiner als dtsmall, so mache keinen neuen Versuch sondern akzeptiere als Not-
c loesung den bestehenden Naeherungswert. Setze dt in diesem Fall als Default
c fuer den kommenden Integrationsschritt auf dtsmall. Setze aber auch das flag
c 'flag_dtsmall', damit gezaehlt werden kann, wie oft dieses Prozedur fuer ein
c bestimmtes Teilchen angewendet werden muss. Ist dies zu oft der Fall, so brich
c diese Trajektorienberechnung ganz ab (-> destiny = code_dtsmall).
c (2. Teil erfolgt weiter unten)
maxErr = Max(maxErr_x,maxErr_v)
if (maxErr.GT.1.) then
if (dt.LT.dtsmall) then ! Fehler immer noch zu gross, obwohl
flag_dtsmall = .true. ! dtsmall schon unterschritten ist
else
!c Bestimme kuerzeren Zeitschritt fuer neuen Versuch (vgl. Referenz):
dt = safety * dt * (maxErr**pShrink)
goto 10
endif
endif
c Nimm die Ergebnisse aus dem dt-Schritt und den beiden dt/2-Schritten und
c berechne damit den neuen Ort und die neue Geschwindigkeit mit Genauigkeit
c fuenfter Ordnung in dt:
x(1) = x(1) + Dx1(1) + xDifferenz(1) / 15.
x(2) = x(2) + Dx1(2) + xDifferenz(2) / 15.
x(3) = x(3) + Dx1(3) + xDifferenz(3) / 15.
v(1) = v(1) + Dv1(1) + vDifferenz(1) / 15.
v(2) = v(2) + Dv1(2) + vDifferenz(2) / 15.
v(3) = v(3) + Dv1(3) + vDifferenz(3) / 15.
c alten Zeitschritt addieren:
t = t + dt
c neuen Zeitschritt so gross wie sinnvoller weise moeglich machen:
if (flag_dtSmall) then
if (n_dtsmall.LT.maxBelowDtSmall) then
dt = dtSmall ! fuer naechsten RK-Schritt
n_dtsmall = n_dtsmall + 1
else
destiny = code_dtsmall ! gib Teilchen verloren
RETURN
endif
else
if (maxErr.GT.errCon) then
dt = safety * dt * (maxErr**pGrow) ! vgl. Referenz
else
dt = 4. * dt ! <- Vergroesserung des Zeitschritts max. um
endif ! Faktor 4!
! pruefen, ob Maximallaenge fuer ersten Testschritt nicht ueberschritten ist:
if (log_confine) dt = min(dt,dl_max/sqrt(v(1)*v(1)+v(2)*v(2)+v(3)*v(3)))
endif
RETURN
c hier folgt der Code fuer 'returnCode_EFeld.NE.0', also fuer den Fall, dass
c bei der Berechnung des Feldes eine Fehlersituation auftrat:
c - Fehler trat auf bei Berechnung des EFeldes am aktuellen Teilchenort oder
c an einem Testort:
999 if (returnCode_EFeld.EQ.2) then
destiny = code_neben_Mappe
elseif (returnCode_EFeld.EQ.3) then
destiny = code_reflektiert
if (v(1).GT.0) then ! in Vorwaertsbewegung -> sollte nicht vorkommen!!
write(*,*)
write(*,*) 'Mappe '//Nr//':'
write(*,*)
write(*,*)' Test-x liegt vor der Mappe!'
write(*,*)' t = ',t
write(*,*)' x = ',x
write(*,*)' v = ',v
write(*,*)' Teilchen-Nr = ',Start_Nr
write(*,*)
endif
elseif (returnCode_EFeld.NE.0) then
write(*,*)
write(*,*) 'SINGLESTEP_RUNGE_KUTTA_'//Nr//': '
write(*,*) 'unallowed value of ''returnCode_EFeld'': ',returnCode_Efeld
write(*,*) '-> STOP'
write(*,*)
STOP
endif
END
c===============================================================================
OPTIONS /EXTEND_SOURCE
SUBROUTINE SINGLESTEP_RUNGE_KUTTA_6(dt,E0,x0,v0, Dx,Dv, *)
c ==========================================================
IMPLICIT NONE
c Diese Subroutine berechnet bei vorgegebenem Zeitschritt einen einzelnen
c Runge-Kutta-Integrationsschritt (4. Ordnung).
c Die Vorgehensweise entspricht den Ausfuehrungen in Kapitel 15.1 der
c NUMERICAL RECIPIES: 'Runge-Kutta Method'.
c Zurueckgegeben werden die errechneten Orts- und Geschwindigkeitsaenderungen
c anstatt direkt der neuen Werte, da sonst vor allem bei den Ortskoordinaten
c Schwierigkeiten auftreten koennen, wenn in der Subroutine 'INTEGRATIONSSTEP_
c RUNGE_KUTTA' aus der Differenz der neuen Werte aus den beiden dt/2- und dem
c dt-Schritt der Fehler abgeschaetzt werden soll (kleine Differenz moeglicher-
c weise grosser Werte).
real Beschl_Faktor
COMMON /Beschl_Faktor/ Beschl_Faktor
real E0(3), x0(3), v0(3) ! Eingangsgroessen
real E1(3), E2(3), E3(3) ! E-Felder an Testorten
real v1(3), v2(3), v3(3) ! Geschwindigkeiten an Testorten
real dt,dt_half,dt_sixth ! zeitl. Aenderung, dt/2, dt/6
real help, help_half, help_sixth ! Hilfsvariable, help/2, help/6
real xTest(3) ! Test-Orte
real Dx(3), Dv(3) ! Ergebnisspeicher
integer i ! Zaehlvariable
c = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = =
dt_half = dt / 2.
dt_sixth = dt / 6.
help = Beschl_Faktor * dt
help_half = help / 2.
help_sixth = help / 6.
do i = 1, 3
xTest(i) = x0(i) + v0(i) * dt_half
v1(i) = v0(i) + E0(i) * help_half
enddo
call EFeld_6(xTest,E1,*999)
do i = 1, 3
xTest(i) = x0(i) + v1(i) * dt_half
v2(i) = v0(i) + E1(i) * help_half
enddo
call EFeld_6(xTest,E2,*999)
do i = 1, 3
xTest(i) = x0(i) + v2(i) * dt
v3(i) = v0(i) + E2(i) * help
enddo
call EFeld_6(xTest,E3,*999)
do i = 1, 3
Dx(i) = (v0(i) + 2.*(v1(i)+v2(i)) + v3(i)) * dt_sixth
Dv(i) = (E0(i) + 2.*(E1(i)+E2(i)) + E3(i)) * help_sixth
enddo
RETURN
999 RETURN 1
END
c===============================================================================
OPTIONS /EXTEND_SOURCE
SUBROUTINE EFeld_6(x,E,*)
c =========================
IMPLICIT NONE
INCLUDE 'accel$sourcedirectory:MAP_DEF_6.INC'
INCLUDE 'accel$sourcedirectory:CALC_FIELD_1.INC'
c - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
c Rechne in Gittereinheiten um:
real_i = (x(1)-xmin) / Dx_
real_j = abs(x(2)) / Dy_
real_k = abs(x(3)) / Dz_
c - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
c Mache die Tests und berechne die Feldstaerke:
c Teste, ob Raumpunkt innerhalb der Potentialmappe liegt:
if (real_j.GT.jmax .OR. real_k.GT.kmax) then
returnCode_EFeld = 2
RETURN 1
elseif (real_i.GT.imax) then
E(1) = 0.
E(2) = 0. ! Spezialbehandlung bei Mappe 6 !
E(3) = 0.
RETURN
elseif (real_i.LT.0.) then
if (real_i.GE.-.1e-6) then
real_i = 0.
else
returnCode_EFeld = 3
RETURN 1
endif
endif
c Bestimme naechstgelegene Stuetzstellen (stuetzstelle_q(n)) und die
c Komponenten des Abstands-Gittervektors zur allernaechsten Stuetzstelle
c (Abstand_q) sowie deren Betraege:
stuetzstelle_i(1) = nint(real_i)
Abstand_i = real_i - stuetzstelle_i(1) ! Abstand zur naeheren Stuetzstelle
Abstand_i_Betrag = abs(Abstand_i)
if (Abstand_i.gt.0.) then
stuetzstelle_i(2) = stuetzstelle_i(1) + 1
elseif (Abstand_i.lt.0.) then
stuetzstelle_i(2) = stuetzstelle_i(1) - 1
else
stuetzstelle_i(2) = stuetzstelle_i(1)
endif
stuetzstelle_j(1) = nint(real_j)
Abstand_j = real_j - stuetzstelle_j(1)
Abstand_j_Betrag = abs(Abstand_j)
if (Abstand_j.gt.0.) then
stuetzstelle_j(2) = stuetzstelle_j(1) + 1
elseif (Abstand_j.lt.0.) then
stuetzstelle_j(2) = stuetzstelle_j(1) - 1
else
stuetzstelle_j(2) = stuetzstelle_j(1)
endif
stuetzstelle_k(1) = nint(real_k)
Abstand_k = real_k - stuetzstelle_k(1)
Abstand_k_Betrag = abs(Abstand_k)
if (Abstand_k.gt.0.) then
stuetzstelle_k(2) = stuetzstelle_k(1) + 1
elseif (Abstand_k.lt.0.) then
stuetzstelle_k(2) = stuetzstelle_k(1) - 1
else
stuetzstelle_k(2) = stuetzstelle_k(1)
endif
c...............................................................................
c Berechnen des elektrischen Feldes:
c ----------------------------------
c
c In dieser Version wird nicht mehr vorausgesetzt, dass das Potential auf den
c Mappenraendern Null ist!
c Bei der Berechnung der Feldstaerke ist angenommen, dass die xy-Ebene (k==0)
c und die xz-Ebene (j==0) Symmetrieebenen sind:
c
c map(i,-j,k) == map(i,j,k).
c map(i,j,-k) == map(i,j,k).
c
c Entlang j=0 ist also E(2)=0, entlang k=0 ist E(3)=0.
c
c (In der vorliegenden Version ist map(i,j,k) durch
c map( k*(jmax+1)*(imax+1) + j*(imax+1) + i) =
c map( (k*(jmax+1) + j)*(imax+1) + i)
c zu ersetzen!)
c (i,j,k laufen von 0 weg, ebenso wie die Indizierung von 'map')
c...............................................................................
c Berechne in den beiden naechstgelegenen k-Ebenen die x-Komponente der Feld-
c staerke. Danach berechne tatsaechlichen Wert aus linearer Interpolation. Um
c die Feldstaerken in den einzelnen k-Ebenen zu bekommen, interpoliere jeweils
c linear zwischen den Werten auf den beiden naechstgelegenen j-Ketten der
c jeweiligen k-Ebene:
i = stuetzstelle_i(1)
do m = 1, 2
k = stuetzstelle_k(m)
do n = 1, 2
j = stuetzstelle_j(n)
ihelp = (k*(jmax+1)+ j)*(imax+1) + i
if (i.EQ.imax) then
c E__(n) = map(imax-1,j,k) - map(imax,j,k)
E__(n) = map(ihelp-1) - map(ihelp )
elseif (i.GT.0) then
c E__(n) = (-0.5+Abstand_i)*(map(i,j,k)-map(i-1,j,k))
c + + ( 0.5+Abstand_i)*(map(i,j,k)-map(i+1,j,k))
E__(n) = (-0.5+Abstand_i)*(map(ihelp)-map(ihelp-1))
+ + ( 0.5+Abstand_i)*(map(ihelp)-map(ihelp+1))
else
c E__(n) = map(0,j,k) - map(1,j,k)
E__(n) = map(ihelp) - map(ihelp+1)
endif
enddo
E_(m) = E__(1) + Abstand_j_Betrag*(E__(2)-E__(1))
enddo
E(1) = E_(1) + Abstand_k_Betrag*(E_(2)-E_(1))
E(1) = E(1) / Dx_ ! Reskalierung entsprechend x-Gitterkonstanten
c Berechne die y-Komponente der Feldstaerke:
j = stuetzstelle_j(1)
do m = 1, 2
k = stuetzstelle_k(m)
do n = 1, 2
i = stuetzstelle_i(n)
ihelp = (k*(jmax+1)+ j)*(imax+1) + i
if (j.EQ.jmax) then
c E__(n) = map(i,jmax-1,k) - map(i,jmax,k)
E__(n) = map(ihelp-(imax+1)) - map(ihelp)
elseif (j.GT.0) then
c E__(n) = (-0.5+Abstand_j)*(map(i,j,k)-map(i,j-1,k))
c + + ( 0.5+Abstand_j)*(map(i,j,k)-map(i,j+1,k))
E__(n) = (-0.5+Abstand_j)*(map(ihelp)-map(ihelp-(imax+1)))
+ + ( 0.5+Abstand_j)*(map(ihelp)-map(ihelp+(imax+1)))
else ! j=0 -> map(i,j-1,k) = map(i,j+1,k) == map(i,1,k)
c E__(n) = 2.0*Abstand_j*(map(i,0,k)-map(i,1,k))
E__(n) = 2.0*Abstand_j*(map(ihelp)-map(ihelp+(imax+1)))
endif
enddo
E_(m) = E__(1) + Abstand_i_Betrag*(E__(2)-E__(1))
enddo
E(2) = E_(1) + Abstand_k_Betrag*(E_(2)-E_(1))
E(2) = E(2) / Dy_ ! Reskalierung entsprechend y-Gitterkonstanten
if (x(2).LT.0) E(2) = -E(2)
c Berechne die z-Komponente der Feldstaerke:
k = stuetzstelle_k(1)
do m = 1, 2
j = stuetzstelle_j(m)
do n = 1, 2
i = stuetzstelle_i(n)
ihelp = (k*(jmax+1)+ j)*(imax+1) + i
if (k.EQ.kmax) then
c E__(n)= map(i,j,kmax-1) - map(i,j,kmax)
E__(n) = map(ihelp-(jmax+1)*(imax+1)) - map(ihelp)
elseif (k.GT.0) then
c E__(n) = (-0.5+Abstand_k)*(map(i,j,k)-map(i,j,k-1))
c + + ( 0.5+Abstand_k)*(map(i,j,k)-map(i,j,k+1))
E__(n) = (-0.5+Abstand_k)*(map(ihelp)-map(ihelp-(jmax+1)*(imax+1)))
+ + ( 0.5+Abstand_k)*(map(ihelp)-map(ihelp+(jmax+1)*(imax+1)))
else ! k=0 -> map(i,j,k-1) = map(i,j,k+1) == map(i,j,1)
c E__(n) = 2.0*Abstand_k*(map(i,j,0)-map(i,j,1))
E__(n) = 2.0*Abstand_k*(map(ihelp)-map(ihelp+(jmax+1)*(imax+1)))
endif
enddo
E_(m) = E__(1) + Abstand_i_Betrag*(E__(2)-E__(1))
enddo
E(3) = E_(1) + Abstand_j_Betrag*(E_(2)-E_(1))
E(3) = E(3) / Dz_ ! Reskalierung entsprechend z-Gitterkonstanten
if (x(3).LT.0) E(3) = -E(3)
cd write(18,*)'x,E = ',x,E
END
c===============================================================================
OPTIONS /EXTEND_SOURCE
SUBROUTINE ERROR_MESSAGE_AND_STOP(text1,text2)
c ==============================================
character*(*) text1,text2
write(*,*)
write(*,*) ' ERROR MESSAGE FROM SUBROUTINE '//text1
write(*,*)
write(*,*) ' >>> '//text2
write(*,*)
write(*,*) ' -> STOP'
write(*,*)
call exit
END
c===============================================================================