module moduleMeshOutputText use moduleOutput, only: prefix, & informFileCreation, & formatFileName, & generateFilePath contains subroutine writeSpeciesOutput(self, fileID, speciesIndex) use moduleMesh use moduleOutput, only: calculateOutput, outputFormat use moduleRefParam, only: L_ref implicit none class(meshParticles), INTENT(in):: self integer, intent(in):: fileID integer, intent(in):: speciesIndex real(8):: r(1:3) type(outputFormat):: output integer:: n do n = 1, self%numNodes r = self%nodes(n)%obj%getCoordinates() call calculateOutput(self%nodes(n)%obj%output(speciesIndex), output, self%nodes(n)%obj%v, species(speciesIndex)%obj) write(fileID, '(5(ES0.6E3,","),ES0.6E3)') r(1)*L_ref, output%density, output%velocity, output%temperature end do end subroutine writeSpeciesOutput subroutine writeCollOutput(self, fileID) use moduleMesh use moduleCollisions implicit none class(meshGeneric), intent(in):: self integer, intent(in):: fileID integer:: n, k, c do n = 1, self%numCells write(fileID, '(I0)', advance='no') n do k = 1, nCollPairs do c = 1, interactionMatrix(k)%amount write(fileID, '(",",I0)', advance='no') self%cells(n)%obj%tallyColl(k)%tally(c) end do end do write(fileID, *) end do end subroutine writeCollOutput subroutine writeEMOutput(self, fileID) use moduleMesh use moduleRefParam, only: L_ref, Volt_ref, B_ref, EF_ref implicit none class(meshParticles), intent(in):: self integer, intent(in):: fileID integer:: n, c real(8):: r(1:3), Xi(1:3) do n = 1, self%numNodes r = self%nodes(n)%obj%getCoordinates() if (n == self%numNodes) then Xi = (/ 1.D0, 0.D0, 0.D0 /) c = self%numNodes - 1 else Xi = (/ 0.D0, 0.D0, 0.D0 /) c = n end if associate(output => self%nodes(n)%obj%emData) write(fileID, '(7(ES0.6E3,","),ES0.6E3)') r(1)*L_ref, & output%phi*Volt_ref, & self%cells(c)%obj%gatherElectricField(Xi)*EF_ref, & output%B*B_ref end associate end do end subroutine writeEMOutput subroutine writeAverage(self, fileIDMean, & fileIDDeviation, & speciesIndex) use moduleMesh use moduleAverage use moduleRefParam, only: L_ref implicit none class(meshParticles), intent(in):: self integer, intent(in):: fileIDMean, fileIDDeviation INTEGER, intent(in):: speciesIndex real(8):: r(1:3) type(outputFormat):: outputMean type(outputFormat):: outputDeviation integer:: n do n = 1, self%numNodes r = self%nodes(n)%obj%getCoordinates() call calculateOutput(averageScheme(n)%mean%output(speciesIndex), outputMean, & self%nodes(n)%obj%v, species(speciesIndex)%obj) write(fileIDMean, '(5(ES0.6E3,","),ES0.6E3)') r(1)*L_ref, outputMean%density, outputMean%velocity, outputMean%temperature call calculateOutput(averageScheme(n)%deviation%output(speciesIndex), outputDeviation, & self%nodes(n)%obj%v, species(speciesIndex)%obj) write(fileIDDeviation, '(5(ES0.6E3,","),ES0.6E3)') r(1)*L_ref, outputDeviation%density, outputDeviation%velocity, outputDeviation%temperature end do end subroutine writeAverage subroutine printOutputText(self) use moduleMesh use moduleSpecies use moduleCaseParam, ONLY: timeStep implicit none class(meshParticles), intent(in):: self INTEGER:: s, fileID character(:), allocatable:: fileName fileID = 60 do s = 1, nSpecies fileName = formatFileName(prefix, species(s)%obj%name, 'csv', timeStep) call informFileCreation(fileName) open (fileID, file = generateFilePath(fileName)) write(fileID, '(5(A,","),A)') '"Position (m)"', & '"Density (m^-3)"', & '"Velocity (m s^-1):0"', 'Velocity (m s^-1):1"', 'Velocity (m s^-1):2"', & '"Temperature (K)"' call writeSpeciesOutput(self, fileID, s) close(fileID) end do end subroutine printOutputText subroutine printCollText(self) use moduleMesh use moduleCaseParam, only: timeStep use moduleOutput, only: collOutput implicit none class(meshGeneric), intent(in):: self integer:: fileID character(:), allocatable:: fileName integer:: k, c character (len=2):: cString fileID = 62 if (collOutput) then fileName = formatFileName(prefix, 'Collisions', 'csv', timeStep) call informFileCreation(fileName) open (fileID, file = generateFilePath(fileName)) write(fileID, '(A)', advance='no') "Cell" do k = 1, nCollPairs do c = 1, interactionMatrix(k)%amount write(cString, "(I2)") c write(fileID, '(",",A)', advance='no') 'Pair ' // interactionMatrix(k)%sp_i%name // '-' // interactionMatrix(k)%sp_j%name // ' collision ' // cString end do end do write(fileID, *) call writeCollOutput(self, fileID) close(fileID) end if end subroutine printCollText subroutine printEMText(self) use moduleMesh use moduleCaseParam, only: timeStep use moduleOutput, only: emOutput implicit none class(meshParticles), intent(in):: self integer:: fileID character(:), allocatable:: fileName fileID = 64 if (emOutput) then fileName = formatFileName(prefix, 'EMField', 'csv', timeStep) call informFileCreation(fileName) open (fileID, file = generateFilePath(fileName)) write(fileID, '(8(A,","),A)') '"Position (m)"', & '"Potential (V)"', & '"Electric Field (V m^-1):0"', 'Electric Field (V m^-1):1"', 'Electric Field (V m^-1):2"', & '"Magnetic Field (T):0"', 'Magnetic Field (T):1"', 'Magnetic Field (T):2"' call writeEMOutput(self, fileID) close(fileID) end if end subroutine printEMText subroutine printAverageText(self) use moduleMesh use moduleSpecies implicit none class(meshParticles), intent(in):: self integer:: s integer:: fileIDMean, fileIDDeviation character(:), allocatable:: fileNameMean, fileNameDeviation fileIDMean = 66 fileIDDeviation = 67 do s = 1, nSpecies fileNameMean = formatFileName('Average_mean', species(s)%obj%name, 'csv') call informFileCreation(fileNameMean) open (fileIDMean, file = generateFilePath(fileNameMean)) write(fileIDMean, '(5(A,","),A)') '"Position (m)"', & '"Density, mean (m^-3)"', & '"Velocity, mean (m s^-1):0"', '"Velocity (m s^-1):1"', '"Velocity (m s^-1):2"', & '"Temperature, mean (K)"' fileNameDeviation = formatFileName('Average_deviation', species(s)%obj%name, 'csv') call informFileCreation(fileNameDeviation) open (fileIDDeviation, file = generateFilePath(fileNameDeviation)) write(fileIDDeviation, '(5(A,","),A)') '"Position (m)"', & '"Density, deviation (m^-3)"', & '"Velocity, deviation (m s^-1):0"', 'Velocity (m s^-1):1"', 'Velocity (m s^-1):2"', & '"Temperature, deviation (K)"' call writeAverage(self, fileIDMean, fileIDDeviation, s) close(fileIDMean) close(fileIDDeviation) end do end subroutine printAverageText end module moduleMeshOutputText