constitutive_nonlocal:

- reworked contribution of immobile dislocation density for rate equations
- flux is now calculated on the basis of interpolated velocities and densities at the interface; both incoming and outgoing fluxes are considered, so every material point only changes his own dotState
- dislocation velocity is now globally defined and calculated by subroutine constitutive_nonlocal_kinetics; the subroutine is called inside _LpAndItsTangent as well as _microstructure; therefore, microstructure now needs Tstar_v as additional input; in the future one should perhaps create a subroutine constitutive_kinetics that calls constitutive_nonlocal_kinetics separately, to clearly distinguish between microstructural and kinetic variables 
- better use flux density vector as output variable instead of scalar flux values for each interface
- added output variables internal and external resolved stress

crystallite:
- added flag to force local stiffness calculation in case of nonlocal model
- misorientation angle is explicitly set to zero when no neighbor can be found

debug:
- added flag "selectiveDebugger" that is used when debugging statements should only affect a specific element, ip and grain; these are specified with the new variables debug_e, debug_i and debug_g
- debugger can now be used in its original sense
This commit is contained in:
Christoph Kords 2010-02-17 13:21:36 +00:00
parent a2037db9c6
commit 97db70cf23
7 changed files with 432 additions and 415 deletions

View File

@ -73,7 +73,10 @@ subroutine CPFEM_general(mode, ffn, ffn1, Temperature, dt, element, IP, cauchySt
use numerics, only: numerics_init, &
relevantStrain, &
iJacoStiffness
use debug, only: debug_init
use debug, only: debug_init, &
debug_g, &
debug_i, &
debug_e
use FEsolving, only: FE_init, &
parallelExecution, &
outdatedFFN1, &
@ -81,6 +84,7 @@ subroutine CPFEM_general(mode, ffn, ffn1, Temperature, dt, element, IP, cauchySt
cycleCounter, &
theInc, &
theTime, &
theDelta, &
FEsolving_execElem, &
FEsolving_execIP
use math, only: math_init, &
@ -205,8 +209,8 @@ subroutine CPFEM_general(mode, ffn, ffn1, Temperature, dt, element, IP, cauchySt
!$OMP CRITICAL (write2out)
write(6,*)
write(6,*) '#####################################'
write(6,'(a10,x,f8.4,x,a10,x,i4,x,a10,x,i3,x,a16,x,i2,x,a16,x,i2)') &
'theTime',theTime,'theInc',theInc,'cycleCounter',cycleCounter,'computationMode',mode
write(6,'(a10,x,f8.4,x,a10,x,f8.4,x,a10,x,i6,x,a10,x,i3,x,a16,x,i2,x,a16,x,i2)') &
'theTime',theTime,'theDelta',theDelta,'theInc',theInc,'cycleCounter',cycleCounter,'computationMode',mode
write(6,*) '#####################################'
call flush (6)
!$OMP END CRITICAL (write2out)
@ -228,7 +232,8 @@ subroutine CPFEM_general(mode, ffn, ffn1, Temperature, dt, element, IP, cauchySt
k = 1:mesh_NcpElems ) &
constitutive_state0(i,j,k)%p = constitutive_state(i,j,k)%p ! microstructure of crystallites
!$OMP CRITICAL (write2out)
write(6,'(a10,/,4(3(e20.8,x),/))') 'aged state',constitutive_state(1,1,1)%p
write(6,'(a,3(x,i4),/,4(3(e20.8,x),/))') 'aged state at', debug_g, debug_i, debug_e, &
constitutive_state(debug_g,debug_i,debug_e)%p
!$OMP END CRITICAL (write2out)
do k = 1,mesh_NcpElems
do j = 1,mesh_maxNips

View File

@ -45,7 +45,7 @@ subroutine constitutive_init()
!* Module initialization *
!**************************************
use prec, only: pReal,pInt
use debug, only: debugger
use debug, only: debugger, selectiveDebugger, debug_e, debug_i, debug_g
use IO, only: IO_error, IO_open_file, IO_open_jobFile
use mesh, only: mesh_maxNips,mesh_NcpElems,mesh_element,FE_Nips
use material
@ -125,7 +125,7 @@ subroutine constitutive_init()
myNgrains = homogenization_Ngrains(mesh_element(3,e))
do i = 1,FE_Nips(mesh_element(2,e)) ! loop over IPs
do g = 1,myNgrains ! loop over grains
debugger = (e == 1 .and. i == 1 .and. g == 1)
selectiveDebugger = (e == debug_e .and. i == debug_i .and. g == debug_g)
myInstance = phase_constitutionInstance(material_phase(g,i,e))
select case(phase_constitution(material_phase(g,i,e)))
@ -269,7 +269,7 @@ return
endfunction
subroutine constitutive_microstructure(Temperature,Fe,Fp,ipc,ip,el)
subroutine constitutive_microstructure(Temperature,Tstar_v,Fe,Fp,ipc,ip,el)
!*********************************************************************
!* This function calculates from state needed variables *
!* INPUT: *
@ -294,6 +294,7 @@ subroutine constitutive_microstructure(Temperature,Fe,Fp,ipc,ip,el)
!* Definition of variables
integer(pInt), intent(in) :: ipc,ip,el
real(pReal), intent(in) :: Temperature
real(pReal), dimension(6) :: Tstar_v
real(pReal), dimension(3,3,homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems), intent(in) :: Fe, Fp
select case (phase_constitution(material_phase(ipc,ip,el)))
@ -308,7 +309,7 @@ real(pReal), dimension(3,3,homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems)
call constitutive_dislotwin_microstructure(Temperature,constitutive_state,ipc,ip,el)
case (constitutive_nonlocal_label)
call constitutive_nonlocal_microstructure(constitutive_state, Temperature, Fe, Fp, ipc, ip, el)
call constitutive_nonlocal_microstructure(constitutive_state, Temperature, Tstar_v, Fe, Fp, ipc, ip, el)
end select

File diff suppressed because it is too large Load Diff

View File

@ -249,6 +249,10 @@ subroutine crystallite_stressAndItsTangent(updateJaco)
nState, &
nCryst
use debug, only: debugger, &
selectiveDebugger, &
debug_e, &
debug_i, &
debug_g, &
debug_CrystalliteLoopDistribution, &
debug_CrystalliteStateLoopDistribution, &
debug_StiffnessStateLoopDistribution
@ -338,6 +342,7 @@ subroutine crystallite_stressAndItsTangent(updateJaco)
logical, dimension(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems) :: &
storedConvergenceFlag
logical, dimension(3,3) :: mask
logical forceLocalStiffnessCalculation = .true. ! flag indicating that stiffness calculation is always done locally
! ------ initialize to starting condition ------
@ -388,9 +393,9 @@ subroutine crystallite_stressAndItsTangent(updateJaco)
myNgrains = homogenization_Ngrains(mesh_element(3,e))
do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e) ! iterate over IPs of this element to be processed
do g = 1,myNgrains
! debugger = (e == 1 .and. i == 1 .and. g == 1)
selectiveDebugger = (e == debug_e .and. i == debug_i .and. g == debug_g)
if (crystallite_converged(g,i,e)) then
if (debugger) then
if (selectiveDebugger) then
!$OMP CRITICAL (write2out)
write(6,'(a21,f10.8,a32,f10.8,a35)') 'winding forward from ', &
crystallite_subFrac(g,i,e),' to current crystallite_subfrac ', &
@ -422,7 +427,7 @@ subroutine crystallite_stressAndItsTangent(updateJaco)
crystallite_Lp(:,:,g,i,e) = crystallite_subLp0(:,:,g,i,e) ! ...plastic velocity grad
constitutive_state(g,i,e)%p = constitutive_subState0(g,i,e)%p ! ...microstructure
crystallite_Tstar_v(:,g,i,e) = crystallite_subTstar0_v(:,g,i,e) ! ...2nd PK stress
if (debugger) then
if (selectiveDebugger) then
!$OMP CRITICAL (write2out)
write(6,'(a78,f10.8)') 'cutback step in crystallite_stressAndItsTangent with new crystallite_subStep: ',&
crystallite_subStep(g,i,e)
@ -463,7 +468,8 @@ subroutine crystallite_stressAndItsTangent(updateJaco)
do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e) ! iterate over IPs of this element to be processed
do g = 1,myNgrains
if (crystallite_todo(g,i,e)) then ! all undone crystallites
call constitutive_microstructure(crystallite_Temperature(g,i,e), crystallite_Fe, crystallite_Fp, g, i, e) ! update dependent state variables to be consistent with basic states
call constitutive_microstructure(crystallite_Temperature(g,i,e), crystallite_Tstar_v(:,g,i,e), crystallite_Fe, &
crystallite_Fp, g, i, e) ! update dependent state variables to be consistent with basic states
constitutive_previousDotState2(g,i,e)%p = 0.0_pReal
constitutive_previousDotState(g,i,e)%p = 0.0_pReal
constitutive_dotState(g,i,e)%p = 0.0_pReal ! zero out dotStates
@ -475,7 +481,7 @@ subroutine crystallite_stressAndItsTangent(updateJaco)
myNgrains = homogenization_Ngrains(mesh_element(3,e))
do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e) ! iterate over IPs of this element to be processed
do g = 1,myNgrains
! debugger = (e == 1 .and. i == 1 .and. g == 1)
! selectiveDebugger = (e == debug_e .and. i == debug_i .and. g == debug_g)
if (crystallite_todo(g,i,e)) then ! all undone crystallites
call constitutive_collectDotState(crystallite_Tstar_v(:,g,i,e), crystallite_subTstar0_v(:,g,i,e), &
crystallite_Fe, crystallite_Fp, crystallite_Temperature(g,i,e), &
@ -491,7 +497,7 @@ subroutine crystallite_stressAndItsTangent(updateJaco)
myNgrains = homogenization_Ngrains(mesh_element(3,e))
do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e) ! iterate over IPs of this element to be processed
do g = 1,myNgrains
! debugger = (e == 1 .and. i == 1 .and. g == 1)
! selectiveDebugger = (e == debug_e .and. i == debug_i .and. g == debug_g)
if (crystallite_todo(g,i,e)) then ! all undone crystallites
crystallite_stateConverged(g,i,e) = crystallite_updateState(g,i,e) ! update state
crystallite_temperatureConverged(g,i,e) = crystallite_updateTemperature(g,i,e) ! update temperature
@ -520,7 +526,7 @@ subroutine crystallite_stressAndItsTangent(updateJaco)
myNgrains = homogenization_Ngrains(mesh_element(3,e))
do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e) ! iterate over IPs of this element to be processed
do g = 1,myNgrains
! debugger = (e == 1 .and. i == 1 .and. g == 1)
! selectiveDebugger = (e == debug_e .and. i == debug_i .and. g == debug_g)
if (crystallite_todo(g,i,e)) & ! all undone crystallites
crystallite_onTrack(g,i,e) = crystallite_integrateStress(g,i,e)
enddo; enddo; enddo
@ -571,7 +577,7 @@ subroutine crystallite_stressAndItsTangent(updateJaco)
myNgrains = homogenization_Ngrains(mesh_element(3,e))
do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e) ! iterate over IPs of this element to be processed
do g = 1,myNgrains
! debugger = (e == 1 .and. i == 1 .and. g == 1)
! selectiveDebugger = (e == debug_e .and. i == debug_i .and. g == debug_g)
if (crystallite_todo(g,i,e)) then ! all undone crystallites
call constitutive_collectDotState(crystallite_Tstar_v(:,g,i,e), crystallite_subTstar0_v(:,g,i,e), &
crystallite_Fe, crystallite_Fp, crystallite_Temperature(g,i,e), &
@ -594,7 +600,7 @@ subroutine crystallite_stressAndItsTangent(updateJaco)
myNgrains = homogenization_Ngrains(mesh_element(3,e))
do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e) ! iterate over IPs of this element to be processed
do g = 1,myNgrains
! debugger = (e == 1 .and. i == 1 .and. g == 1)
! selectiveDebugger = (e == debug_e .and. i == debug_i .and. g == debug_g)
if (crystallite_todo(g,i,e)) then ! all undone crystallites
crystallite_stateConverged(g,i,e) = crystallite_updateState(g,i,e) ! update state
crystallite_temperatureConverged(g,i,e) = crystallite_updateTemperature(g,i,e) ! update temperature
@ -615,7 +621,11 @@ subroutine crystallite_stressAndItsTangent(updateJaco)
!$OMP CRITICAL (write2out)
write(6,*) count(crystallite_converged(:,:,:)),'grains converged after state integration no.', NiterationState
write(6,*)
!write(6,'(8(L,x))') crystallite_converged(:,:,:)
! write(6,'(8(L,x))') crystallite_converged(:,:,:)
! do e = FEsolving_execElem(1),FEsolving_execElem(2)
! if (any(.not. crystallite_converged(:,:,e))) &
! write(6,'(i4,8(x,L))') e, crystallite_converged(:,:,e)
! enddo
!$OMPEND CRITICAL (write2out)
endif
if (any(.not. crystallite_converged .and. .not. crystallite_localConstitution)) & ! any non-local not yet converged?
@ -684,24 +694,24 @@ subroutine crystallite_stressAndItsTangent(updateJaco)
storedP = crystallite_P
storedConvergenceFlag = crystallite_converged
if (all(crystallite_localConstitution) .or. theInc < 2) then ! all grains have local constitution, so local convergence of perturbed grain is sufficient
if (all(crystallite_localConstitution) .or. theInc < 2 .or. forceLocalStiffnessCalculation) then ! all grains have local constitution, so local convergence of perturbed grain is sufficient
!$OMP PARALLEL DO
do e = FEsolving_execElem(1),FEsolving_execElem(2) ! iterate over elements to be processed
myNgrains = homogenization_Ngrains(mesh_element(3,e))
do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e) ! iterate over IPs of this element to be processed
do g = 1,myNgrains
! debugger = (e==1 .and. i==1 .and. g==1)
! selectiveDebugger = (e == debug_e .and. i == debug_i .and. g == debug_g)
if (crystallite_requested(g,i,e)) then ! first check whether is requested at all!
if (crystallite_converged(g,i,e)) then ! grain converged in above iteration
if (debugger) then
if (selectiveDebugger) then
!$OMP CRITICAL (write2out)
write (6,*) '#############'
write (6,*) 'central solution of cryst_StressAndTangent'
write (6,*) '#############'
write (6,'(a,/,3(3(f12.4,x)/))') ' P of 1 1 1',storedP(1:3,:,g,i,e)/1e6
write (6,'(a,/,3(3(f12.8,x)/))') ' Fp of 1 1 1',storedFp(1:3,:,g,i,e)
write (6,'(a,/,3(3(f12.8,x)/))') ' Lp of 1 1 1',storedLp(1:3,:,g,i,e)
write (6,'(a8,3(x,i4),/,3(3(f12.4,x)/))') ' P of', g, i, e, storedP(1:3,:,g,i,e)/1e6
write (6,'(a8,3(x,i4),/,3(3(f12.8,x)/))') ' Fp of', g, i, e, storedFp(1:3,:,g,i,e)
write (6,'(a8,3(x,i4),/,3(3(f12.8,x)/))') ' Lp of', g, i, e, storedLp(1:3,:,g,i,e)
!$OMPEND CRITICAL (write2out)
endif
@ -711,10 +721,10 @@ subroutine crystallite_stressAndItsTangent(updateJaco)
do k = 1,3 ! perturbation...
do l = 1,3 ! ...components to the positive direction
crystallite_subF(k,l,g,i,e) = crystallite_subF(k,l,g,i,e) + myPert ! perturb single component (either forward or backward)
if (debugger) then
if (selectiveDebugger) then
!$OMP CRITICAL (write2out)
write (6,'(i1,x,i1)') k,l
write (6,'(a,/,3(3(f12.6,x)/))') 'pertF of 1 1 1',crystallite_subF(1:3,:,g,i,e)
write (6,'(a8,3(x,i4),/,3(3(f12.6,x)/))') 'pertF of', g, i, e, crystallite_subF(1:3,:,g,i,e)
!$OMPEND CRITICAL (write2out)
endif
onTrack = .true.
@ -733,12 +743,13 @@ subroutine crystallite_stressAndItsTangent(updateJaco)
temperatureConverged = crystallite_updateTemperature(g,i,e) ! update temperature
converged = stateConverged .and. temperatureConverged
endif
if (debugger) then
if (selectiveDebugger) then
!$OMP CRITICAL (write2out)
write (6,*) '-------------'
write (6,'(a,x,l,x,l)') 'ontrack + converged:',onTrack,converged
write (6,'(a,/,3(3(f12.4,x)/))') 'pertP/MPa of 1 1 1',crystallite_P(1:3,:,g,i,e)/1e6
write (6,'(a,/,3(3(f12.4,x)/))') 'DP/MPa of 1 1 1',(crystallite_P(1:3,:,g,i,e)-storedP(1:3,:,g,i,e))/1e6
write (6,'(a12,3(x,i4),/,3(3(f12.4,x)/))') 'pertP/MPa of', g, i, e, crystallite_P(1:3,:,g,i,e)/1e6
write (6,'(a12,3(x,i4),/,3(3(f12.4,x)/))') 'DP/MPa of', g, i, e, &
(crystallite_P(1:3,:,g,i,e)-storedP(1:3,:,g,i,e))/1e6
!$OMPEND CRITICAL (write2out)
endif
enddo
@ -918,7 +929,8 @@ endsubroutine
constitutive_state, &
constitutive_relevantState, &
constitutive_microstructure
use debug, only: debugger
use debug, only: debugger, &
selectiveDebugger
use FEsolving, only: cycleCounter, theInc
!*** input variables ***!
@ -957,13 +969,14 @@ endsubroutine
! update the microstructure
constitutive_state(g,i,e)%p(1:mySize) = constitutive_state(g,i,e)%p(1:mySize) - residuum
call constitutive_microstructure(crystallite_Temperature(g,i,e), crystallite_Fe, crystallite_Fp, g, i, e)
call constitutive_microstructure(crystallite_Temperature(g,i,e), crystallite_Tstar_v(:,g,i,e), crystallite_Fe, crystallite_Fp, &
g, i, e)
! setting flag to true if state is below relative tolerance, otherwise set it to false
crystallite_updateState = all( constitutive_state(g,i,e)%p(1:mySize) < constitutive_relevantState(g,i,e)%p(1:mySize) &
.or. abs(residuum) < rTol_crystalliteState*abs(constitutive_state(g,i,e)%p(1:mySize)))
if (debugger) then
if (selectiveDebugger) then
!$OMP CRITICAL (write2out)
if (crystallite_updateState) then
write(6,*) '::: updateState converged',g,i,e
@ -1073,6 +1086,7 @@ endsubroutine
iJacoLpresiduum, &
relevantStrain
use debug, only: debugger, &
selectiveDebugger, &
debug_cumLpCalls, &
debug_cumLpTicks, &
debug_StressLoopDistribution
@ -1214,13 +1228,12 @@ LpLoop: do
debug_cumLpCalls = debug_cumLpCalls + 1_pInt
debug_cumLpTicks = debug_cumLpTicks + tock-tick
if (tock < tick) debug_cumLpTicks = debug_cumLpTicks + maxticks
if (debugger) then
if (selectiveDebugger) then
!$OMP CRITICAL (write2out)
write(6,*) '::: integrateStress at ' ,g,i,e, ' ; iteration ', NiterationStress
write(6,*)
write(6,'(a,/,3(3(f20.7,x)/))') 'Lp_constitutive', Lp_constitutive
write(6,'(a,/,3(3(f20.7,x)/))') 'Lpguess', Lpguess
! call flush(6)
!$OMPEND CRITICAL (write2out)
endif
@ -1332,7 +1345,7 @@ LpLoop: do
! set return flag to true
crystallite_integrateStress = .true.
if (debugger) then
if (selectiveDebugger) then
!$OMP CRITICAL (write2out)
write(6,*) '::: integrateStress converged at ',g,i,e,' ; iteration ', NiterationStress
write(6,*)
@ -1468,11 +1481,12 @@ logical error
crystallite_R(:,:,1,neighboring_i,neighboring_e), &
symmetryType) ! calculate misorientation
else
crystallite_misorientation(4,n,1,i,e) = 400.0_pReal
! call IO_warning(-1, e, i, g, 'couldnt calculate misorientation because of two different lattice structures')
else ! for neighbor with different phase
crystallite_misorientation(4,n,1,i,e) = 400.0_pReal ! set misorientation angle to 400
endif
else ! no existing neighbor
crystallite_misorientation(4,n,1,i,e) = 0.0_pReal ! set misorientation angle to zero
endif
enddo
endif

View File

@ -17,6 +17,10 @@
integer(pInt) :: debug_cumLpCalls = 0_pInt
integer(pInt) :: debug_cumDotStateCalls = 0_pInt
integer(pInt) :: debug_cumDotTemperatureCalls = 0_pInt
integer(pInt) :: debug_e = 1_pInt
integer(pInt) :: debug_i = 1_pInt
integer(pInt) :: debug_g = 1_pInt
logical :: selectiveDebugger = .false.
logical :: debugger = .false.
logical :: distribution_init = .false.

View File

@ -220,6 +220,9 @@ subroutine materialpoint_stressAndItsTangent(&
crystallite_stressAndItsTangent, &
crystallite_orientations
use debug, only: debugger, &
selectiveDebugger, &
debug_e, &
debug_i, &
debug_MaterialpointLoopDistribution, &
debug_MaterialpointStateLoopDistribution
@ -278,11 +281,11 @@ subroutine materialpoint_stressAndItsTangent(&
myNgrains = homogenization_Ngrains(mesh_element(3,e))
do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e) ! iterate over IPs of this element to be processed
debugger = (e == 1 .and. i == 1)
selectiveDebugger = (e == debug_e .and. i == debug_i)
! if our materialpoint converged or consists of only one single grain then we are either finished or have to wind forward
if ( materialpoint_converged(i,e) .or. (myNgrains == 1_pInt .and. materialpoint_subStep(i,e) <= 1.0_pReal) ) then
if (debugger) then
if (selectiveDebugger) then
!$OMP CRITICAL (write2out)
write(6,'(a21,f10.8,a34,f10.8,a37,/)') 'winding forward from ', &
materialpoint_subFrac(i,e), ' to current materialpoint_subFrac ', &
@ -321,7 +324,7 @@ subroutine materialpoint_stressAndItsTangent(&
materialpoint_subStep(i,e) = subStepSizeHomog * materialpoint_subStep(i,e) ! crystallite had severe trouble, so do a significant cutback
! <<modified to add more flexibility in cutback>>
if (debugger) then
if (selectiveDebugger) then
!$OMP CRITICAL (write2out)
write(6,'(a82,f10.8,/)') 'cutback step in materialpoint_stressAndItsTangent with new materialpoint_subStep: ',&
materialpoint_subStep(i,e)

View File

@ -149,6 +149,8 @@ constitution nonlocal
(output) delta_dip
(output) shearrate
(output) resolvedstress
(output) resolvedstress_internal
(output) resolvedstress_external
(output) resistance
(output) rho_dot
(output) rho_dot_sgl
@ -160,35 +162,18 @@ constitution nonlocal
(output) rho_dot_dip2sgl
(output) rho_dot_ann_ath
(output) rho_dot_ann_the
(output) rho_dot_flux
(output) rho_dot_flux_ep
(output) rho_dot_flux_en
(output) rho_dot_flux_sp
(output) rho_dot_flux_sn
(output) rho_dot_flux_n1_ep
(output) rho_dot_flux_n1_en
(output) rho_dot_flux_n1_sp
(output) rho_dot_flux_n1_sn
(output) rho_dot_flux_n2_ep
(output) rho_dot_flux_n2_en
(output) rho_dot_flux_n2_sp
(output) rho_dot_flux_n2_sn
(output) rho_dot_flux_n3_ep
(output) rho_dot_flux_n3_en
(output) rho_dot_flux_n3_sp
(output) rho_dot_flux_n3_sn
(output) rho_dot_flux_n4_ep
(output) rho_dot_flux_n4_en
(output) rho_dot_flux_n4_sp
(output) rho_dot_flux_n4_sn
(output) rho_dot_flux_n5_ep
(output) rho_dot_flux_n5_en
(output) rho_dot_flux_n5_sp
(output) rho_dot_flux_n5_sn
(output) rho_dot_flux_n6_ep
(output) rho_dot_flux_n6_en
(output) rho_dot_flux_n6_sp
(output) rho_dot_flux_n6_sn
(output) fluxDensity_edge_pos_x
(output) fluxDensity_edge_pos_y
(output) fluxDensity_edge_pos_z
(output) fluxDensity_edge_neg_x
(output) fluxDensity_edge_neg_y
(output) fluxDensity_edge_neg_z
(output) fluxDensity_screw_pos_x
(output) fluxDensity_screw_pos_y
(output) fluxDensity_screw_pos_z
(output) fluxDensity_screw_neg_x
(output) fluxDensity_screw_neg_y
(output) fluxDensity_screw_neg_z
(output) d_upper_edge
(output) d_upper_screw
(output) d_upper_dot_edge