!!****************************************************************************** !! !! This file is part of the AMUN source code, a program to perform !! Newtonian or relativistic magnetohydrodynamical simulations on uniform or !! adaptive mesh. !! !! Copyright (C) 2008-2021 Grzegorz Kowal !! !! This program is free software: you can redistribute it and/or modify !! it under the terms of the GNU General Public License as published by !! the Free Software Foundation, either version 3 of the License, or !! (at your option) any later version. !! !! This program is distributed in the hope that it will be useful, !! but WITHOUT ANY WARRANTY; without even the implied warranty of !! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the !! GNU General Public License for more details. !! !! You should have received a copy of the GNU General Public License !! along with this program. If not, see . !! !!****************************************************************************** !! !! module: INTEGRALS !! !! This module provides subroutines to calculate and store integrals of !! the conserved variables, and other useful statistics in the integrals file. !! !!****************************************************************************** ! module integrals #ifdef PROFILE ! import external subroutines ! use timers, only : set_timer, start_timer, stop_timer #endif /* PROFILE */ ! module variables are not implicit by default ! implicit none #ifdef PROFILE ! timer indices ! integer, save :: imi, ims #endif /* PROFILE */ ! MODULE PARAMETERS: ! ================= ! ! funit - a file handler to the integrals file; ! sunit - a file handler to the statistics file; ! eunit - a file handler to the errors file; ! iintd - the number of steps between subsequent intervals storing; ! integer(kind=4), save :: funit = 11 integer(kind=4), save :: sunit = 12 integer(kind=4), save :: eunit = 13 integer(kind=4), save :: iintd = 1 ! flag indicating if the files are actually written to disk ! logical, save :: stored = .false. ! the format string ! character(len=32), save :: efmt ! by default everything is private ! private ! declare public subroutines ! public :: initialize_integrals, finalize_integrals public :: store_integrals !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ! contains ! !=============================================================================== !! !!*** PUBLIC SUBROUTINES ***************************************************** !! !=============================================================================== ! !=============================================================================== ! ! subroutine INITIALIZE_INTEGRALS: ! ------------------------------- ! ! Subroutine initializes module INTEGRALS by preparing the integrals file. ! ! Arguments: ! ! store - flag indicating if the files should be actually written to disk; ! irun - job execution counter; ! status - return flag of the procedure execution status; ! !=============================================================================== ! subroutine initialize_integrals(store, irun, status) ! import external variables and subroutines ! use equations , only : pvars, nf use evolution , only : error_control use parameters, only : get_parameter ! local variables are not implicit by default ! implicit none ! subroutine arguments ! logical, intent(in) :: store integer, intent(inout) :: irun integer, intent(out) :: status ! local variables ! character(len=32) :: fname, append, stmp logical :: flag ! !------------------------------------------------------------------------------- ! #ifdef PROFILE ! set timer descriptions ! call set_timer('integrals:: initialization' , imi) call set_timer('integrals:: integrals storing', ims) ! start accounting time for module initialization/finalization ! call start_timer(imi) #endif /* PROFILE */ ! reset the status flag ! status = 0 ! set stored flag ! stored = store ! get the integrals storing interval ! call get_parameter("integrals_interval", iintd) ! make sure storage interval is larger than zero ! iintd = max(1, iintd) ! handle the integral file if store is set ! if (stored) then ! depending on the append parameter, choose the right file initialization for ! the integrals file ! append = "off" call get_parameter("integrals_append", append) select case(trim(append)) case ("on", "ON", "t", "T", "y", "Y", "true", "TRUE", "yes", "YES") write(fname, "('integrals.dat')") inquire(file=fname, exist=flag) case default write(fname, "('integrals_',i2.2,'.dat')") irun flag = .false. end select ! check if the file exists; if not, create a new one, otherwise move to the end ! if (flag .and. irun > 1) then #ifdef __INTEL_COMPILER open(newunit=funit, file=fname, form='formatted', status='old', & position='append', buffered='yes') #else /* __INTEL_COMPILER */ open(newunit=funit, file=fname, form='formatted', status='old', & position='append') #endif /* __INTEL_COMPILER */ write(funit,"('#')") else #ifdef __INTEL_COMPILER open(newunit=funit, file=fname, form='formatted', & status='replace', buffered='yes') #else /* __INTEL_COMPILER */ open(newunit=funit, file=fname, form='formatted', status='replace') #endif /* __INTEL_COMPILER */ end if ! write the integral file header ! write(funit,"('#',a8,13(1x,a18))") 'step', 'time', 'dt' & , 'mass', 'momx', 'momy', 'momz' & , 'ener', 'ekin', 'emag', 'eint' & , 'einj', 'erat', 'arms' write(funit,"('#')") ! depending on the append parameter, choose the right file initialization for ! the statistics file ! append = "off" call get_parameter("statistics_append", append) select case(trim(append)) case ("on", "ON", "t", "T", "y", "Y", "true", "TRUE", "yes", "YES") write(fname, "('statistics.dat')") inquire(file=fname, exist=flag) case default write(fname, "('statistics_',i2.2,'.dat')") irun flag = .false. end select ! check if the file exists; if not, create a new one, otherwise move to the end ! if (flag .and. irun > 1) then #ifdef __INTEL_COMPILER open(newunit=sunit, file=fname, form='formatted', status='old', & position='append', buffered='yes') #else /* __INTEL_COMPILER */ open(newunit=sunit, file=fname, form='formatted', status='old', & position='append') #endif /* __INTEL_COMPILER */ write(sunit,"('#')") else #ifdef __INTEL_COMPILER open(newunit=sunit, file=fname, form='formatted', & status='replace', buffered='yes') #else /* __INTEL_COMPILER */ open(newunit=sunit, file=fname, form='formatted', status='replace') #endif /* __INTEL_COMPILER */ end if ! write the integral file header ! write(sunit,"('#',a8,23(1x,a18))") 'step', 'time' & , 'mean(dens)', 'min(dens)', 'max(dens)' & , 'mean(pres)', 'min(pres)', 'max(pres)' & , 'mean(velo)', 'min(velo)', 'max(velo)' & , 'mean(magn)', 'min(magn)', 'max(magn)' & , 'mean(bpsi)', 'min(bpsi)', 'max(bpsi)' & , 'mean(Mson)', 'min(Mson)', 'max(Mson)' & , 'mean(Malf)', 'min(Malf)', 'max(Malf)' write(sunit,"('#')") ! depending on the append parameter, choose the right file initialization for ! the file to store integration errors ! if (error_control) then append = "off" call get_parameter("errors_append", append) select case(trim(append)) case ("on", "ON", "t", "T", "y", "Y", "true", "TRUE", "yes", "YES") write(fname, "('errors.dat')") inquire(file=fname, exist=flag) case default write(fname, "('errors_',i2.2,'.dat')") irun flag = .false. end select ! check if the file exists; if not, create a new one, otherwise move to the end ! if (flag .and. irun > 1) then #ifdef __INTEL_COMPILER open(newunit=eunit, file=fname, form='formatted', status='old', & position='append', buffered='yes') #else /* __INTEL_COMPILER */ open(newunit=eunit, file=fname, form='formatted', status='old', & position='append') #endif /* __INTEL_COMPILER */ write(funit,"('#')") else #ifdef __INTEL_COMPILER open(newunit=eunit, file=fname, form='formatted', & status='replace', buffered='yes') #else /* __INTEL_COMPILER */ open(newunit=eunit, file=fname, form='formatted', status='replace') #endif /* __INTEL_COMPILER */ end if ! write the integral file header ! write(stmp,"(i9)") nf + 4 efmt = "('#',a8," // trim(adjustl(stmp)) // "(1x,a18))" write(eunit,efmt) 'step', 'time', 'dt_cfl', 'dt_err', & 'max_err', pvars(1:nf) write(eunit,"('#')") ! prepare format for errors ! efmt = "(i9," // trim(adjustl(stmp)) // "(1x,1es18.8e3))" end if ! error_control end if ! store #ifdef PROFILE ! stop accounting time for module initialization/finalization ! call stop_timer(imi) #endif /* PROFILE */ !------------------------------------------------------------------------------- ! end subroutine initialize_integrals ! !=============================================================================== ! ! subroutine FINALIZE_INTEGRALS: ! ----------------------------- ! ! Subroutine finalizes module INTEGRALS by closing the integrals file. ! ! Arguments: ! ! status - an integer flag for error return value; ! !=============================================================================== ! subroutine finalize_integrals(status) use evolution, only : error_control implicit none integer, intent(out) :: status !------------------------------------------------------------------------------- ! #ifdef PROFILE call start_timer(imi) #endif /* PROFILE */ status = 0 ! close the integrals and statistics files ! if (stored) then close(funit) close(sunit) if (error_control) close(eunit) end if #ifdef PROFILE call stop_timer(imi) #endif /* PROFILE */ !------------------------------------------------------------------------------- ! end subroutine finalize_integrals ! !=============================================================================== ! ! subroutine STORE_INTEGRALS: ! -------------------------- ! ! Subroutine calculates the integrals, collects from all processes and ! stores them in the integrals file. ! ! !=============================================================================== ! subroutine store_integrals() use blocks , only : block_meta, block_data, list_data use coordinates , only : ni => ncells, nb, ne use coordinates , only : advol, voli use equations , only : idn, ipr, ivx, ivz, ibx, iby, ibz, ibp use equations , only : ien, imx, imy, imz use equations , only : magnetized, adiabatic_index, csnd use equations , only : errors use evolution , only : error_control use evolution , only : step, time, dt, dth, dte use forcing , only : einj, rinj, arms use helpers , only : flush_and_sync use iso_fortran_env, only : error_unit #ifdef MPI use mpitools , only : reduce_minimum, reduce_maximum, reduce_sum #endif /* MPI */ use workspace , only : work, nwork, work_in_use implicit none logical, save :: first = .true. integer :: n, l, u, nk, status real(kind=8) :: dvol, dvolh type(block_data), pointer :: pdata real(kind=8), dimension(:,:,:), pointer, save :: vel, mag, sqd, tmp integer , parameter :: narr = 16 real(kind=8), dimension(narr) :: inarr, avarr, mnarr, mxarr real(kind=8) , parameter :: eps = epsilon(1.0d+00) real(kind=8) , parameter :: big = huge(1.0d+00) character(len=*), parameter :: loc = 'INTEGRALS:store_integrals()' !------------------------------------------------------------------------------- ! ! return if the storage interval was not reached ! if (mod(step, iintd) > 0) return #ifdef PROFILE call start_timer(ims) #endif /* PROFILE */ if (first) then status = 0 if (work_in_use) then write(error_unit,"('[',a,']: ',a,3i4,a)") trim(loc), & "Workspace is already occupied!" return else n = 4 * ni**NDIMS if (n > nwork) then write(error_unit,"('[',a,']: ',a,3i4,a)") trim(loc), & "Workspace has too be increased!" nwork = n deallocate(work) allocate(work(nwork), stat=status) if (status /= 0) then write(error_unit,"('[',a,']: ',a)") trim(loc), & "Cannot increase the workspace's size!" return end if end if #if NDIMS == 3 nk = ni #else /* NDIMS == 3 */ nk = 1 #endif /* NDIMS == 3 */ n = ni**NDIMS l = 1 u = n vel(1:ni,1:ni,1:nk) => work(l:u) l = l + n u = u + n mag(1:ni,1:ni,1:nk) => work(l:u) l = l + n u = u + n sqd(1:ni,1:ni,1:nk) => work(l:u) l = l + n u = u + n tmp(1:ni,1:ni,1:nk) => work(l:u) end if first = .false. end if ! reset the integrals array ! inarr(:) = 0.0d+00 avarr(:) = 0.0d+00 mnarr(:) = big mxarr(:) = - big ! reset some statistics if they are not used ! if (ipr < 1) then mnarr(2) = 0.0d+00 mxarr(2) = 0.0d+00 end if if (.not. magnetized) then mnarr(4) = 0.0d+00 mxarr(4) = 0.0d+00 mnarr(5) = 0.0d+00 mxarr(5) = 0.0d+00 mnarr(7) = 0.0d+00 mxarr(7) = 0.0d+00 end if if (work_in_use) then write(error_unit,"('[',a,']: ',a,3i4,a)") trim(loc), & "Workspace is already occupied! Corruptions can occur!" end if work_in_use = .true. ! associate the pointer with the first block on the data block list ! pdata => list_data ! iterate over all data blocks on the list ! do while(associated(pdata)) ! obtain the volume elements for the current block ! dvol = advol(pdata%meta%level) dvolh = 0.5d+00 * dvol ! sum up density and momenta components ! #if NDIMS == 3 inarr(1) = inarr(1) + sum(pdata%u(idn,nb:ne,nb:ne,nb:ne)) * dvol inarr(2) = inarr(2) + sum(pdata%u(imx,nb:ne,nb:ne,nb:ne)) * dvol inarr(3) = inarr(3) + sum(pdata%u(imy,nb:ne,nb:ne,nb:ne)) * dvol inarr(4) = inarr(4) + sum(pdata%u(imz,nb:ne,nb:ne,nb:ne)) * dvol #else /* NDIMS == 3 */ inarr(1) = inarr(1) + sum(pdata%u(idn,nb:ne,nb:ne, : )) * dvol inarr(2) = inarr(2) + sum(pdata%u(imx,nb:ne,nb:ne, : )) * dvol inarr(3) = inarr(3) + sum(pdata%u(imy,nb:ne,nb:ne, : )) * dvol inarr(4) = inarr(4) + sum(pdata%u(imz,nb:ne,nb:ne, : )) * dvol #endif /* NDIMS == 3 */ ! sum up total energy ! if (ien > 0) then #if NDIMS == 3 inarr(5) = inarr(5) + sum(pdata%u(ien,nb:ne,nb:ne,nb:ne)) * dvol #else /* NDIMS == 3 */ inarr(5) = inarr(5) + sum(pdata%u(ien,nb:ne,nb:ne, : )) * dvol #endif /* NDIMS == 3 */ end if ! sum up kinetic energy ! #if NDIMS == 3 inarr(6) = inarr(6) + sum((pdata%u(imx,nb:ne,nb:ne,nb:ne)**2 & + pdata%u(imy,nb:ne,nb:ne,nb:ne)**2 & + pdata%u(imz,nb:ne,nb:ne,nb:ne)**2) & / pdata%u(idn,nb:ne,nb:ne,nb:ne)) * dvolh #else /* NDIMS == 3 */ inarr(6) = inarr(6) + sum((pdata%u(imx,nb:ne,nb:ne, : )**2 & + pdata%u(imy,nb:ne,nb:ne, : )**2 & + pdata%u(imz,nb:ne,nb:ne, : )**2) & / pdata%u(idn,nb:ne,nb:ne, : )) * dvolh #endif /* NDIMS == 3 */ ! sum up magnetic energy ! if (magnetized) then #if NDIMS == 3 inarr(7) = inarr(7) + sum(pdata%u(ibx,nb:ne,nb:ne,nb:ne)**2 & + pdata%u(iby,nb:ne,nb:ne,nb:ne)**2 & + pdata%u(ibz,nb:ne,nb:ne,nb:ne)**2) * dvolh #else /* NDIMS == 3 */ inarr(7) = inarr(7) + sum(pdata%u(ibx,nb:ne,nb:ne, : )**2 & + pdata%u(iby,nb:ne,nb:ne, : )**2 & + pdata%u(ibz,nb:ne,nb:ne, : )**2) * dvolh #endif /* NDIMS == 3 */ end if ! sum up the injected energy and injection rate ! inarr( 9) = einj inarr(10) = rinj ! get average, minimum and maximum values of density ! #if NDIMS == 3 tmp(:,:,:) = pdata%q(idn,nb:ne,nb:ne,nb:ne) #else /* NDIMS == 3 */ tmp(:,:,:) = pdata%q(idn,nb:ne,nb:ne, : ) #endif /* NDIMS == 3 */ avarr(1) = avarr(1) + sum(tmp(:,:,:)) * dvol mnarr(1) = min(mnarr(1), minval(tmp(:,:,:))) mxarr(1) = max(mxarr(1), maxval(tmp(:,:,:))) ! get average, minimum and maximum values of pressure ! if (ipr > 0) then #if NDIMS == 3 tmp(:,:,:) = pdata%q(ipr,nb:ne,nb:ne,nb:ne) #else /* NDIMS == 3 */ tmp(:,:,:) = pdata%q(ipr,nb:ne,nb:ne, : ) #endif /* NDIMS == 3 */ avarr(2) = avarr(2) + sum(tmp(:,:,:)) * dvol mnarr(2) = min(mnarr(2), minval(tmp(:,:,:))) mxarr(2) = max(mxarr(2), maxval(tmp(:,:,:))) end if ! get average, minimum and maximum values of velocity amplitude ! #if NDIMS == 3 vel(:,:,:) = sqrt(sum(pdata%q(ivx:ivz,nb:ne,nb:ne,nb:ne)**2, 1)) #else /* NDIMS == 3 */ vel(:,:,:) = sqrt(sum(pdata%q(ivx:ivz,nb:ne,nb:ne, : )**2, 1)) #endif /* NDIMS == 3 */ avarr(3) = avarr(3) + sum(vel(:,:,:)) * dvol mnarr(3) = min(mnarr(3), minval(vel(:,:,:))) mxarr(3) = max(mxarr(3), maxval(vel(:,:,:))) ! get average, minimum and maximum values of magnetic field amplitude, and ! divergence potential ! if (magnetized) then #if NDIMS == 3 mag(:,:,:) = sqrt(sum(pdata%q(ibx:ibz,nb:ne,nb:ne,nb:ne)**2, 1)) #else /* NDIMS == 3 */ mag(:,:,:) = sqrt(sum(pdata%q(ibx:ibz,nb:ne,nb:ne, : )**2, 1)) #endif /* NDIMS == 3 */ avarr(4) = avarr(4) + sum(mag(:,:,:)) * dvol mnarr(4) = min(mnarr(4), minval(mag(:,:,:))) mxarr(4) = max(mxarr(4), maxval(mag(:,:,:))) #if NDIMS == 3 tmp(:,:,:) = pdata%q(ibp,nb:ne,nb:ne,nb:ne) #else /* NDIMS == 3 */ tmp(:,:,:) = pdata%q(ibp,nb:ne,nb:ne, : ) #endif /* NDIMS == 3 */ avarr(5) = avarr(5) + sum(tmp(:,:,:)) * dvol mnarr(5) = min(mnarr(5), minval(tmp(:,:,:))) mxarr(5) = max(mxarr(5), maxval(tmp(:,:,:))) end if ! calculate square root of density ! #if NDIMS == 3 sqd(:,:,:) = sqrt(pdata%q(idn,nb:ne,nb:ne,nb:ne)) #else /* NDIMS == 3 */ sqd(:,:,:) = sqrt(pdata%q(idn,nb:ne,nb:ne, : )) #endif /* NDIMS == 3 */ ! get average, minimum and maximum values of sonic Mach number ! if (ipr > 0) then tmp(:,:,:) = sqd(:,:,:) * vel(:,:,:) / & #if NDIMS == 3 sqrt(adiabatic_index * pdata%q(ipr,nb:ne,nb:ne,nb:ne)) #else /* NDIMS == 3 */ sqrt(adiabatic_index * pdata%q(ipr,nb:ne,nb:ne, : )) #endif /* NDIMS == 3 */ else tmp(:,:,:) = vel(:,:,:) / csnd end if avarr(6) = avarr(6) + sum(tmp(:,:,:)) * dvol mnarr(6) = min(mnarr(6), minval(tmp(:,:,:))) mxarr(6) = max(mxarr(6), maxval(tmp(:,:,:))) ! get average, minimum and maximum values of Alfvénic Mach number ! if (magnetized) then tmp(:,:,:) = sqd(:,:,:) * vel(:,:,:) / max(eps, mag(:,:,:)) avarr(7) = avarr(7) + sum(tmp(:,:,:)) * dvol mnarr(7) = min(mnarr(7), minval(tmp(:,:,:))) mxarr(7) = max(mxarr(7), maxval(tmp(:,:,:))) end if ! associate the pointer with the next block on the list ! pdata => pdata%next end do ! data blocks work_in_use = .false. #ifdef MPI ! sum the integral array from all processes ! call reduce_sum(inarr(1:narr)) ! reduce average, minimum and maximum values ! call reduce_sum(avarr(1:narr)) call reduce_minimum(mnarr(1:narr)) call reduce_maximum(mxarr(1:narr)) #endif /* MPI */ ! calculate the internal energy ! if (ien > 0) inarr(8) = inarr(5) - inarr(6) - inarr(7) ! normalize the averages by the volume of domain ! avarr(:) = avarr(:) * voli ! write down the integrals and statistics to appropriate files ! if (stored) then write(funit,"(i9,13(1x,1es18.8e3))") step, time, dt, inarr(1:10), arms write(sunit,"(i9,23(1x,1es18.8e3))") step, time & , avarr(1), mnarr(1), mxarr(1) & , avarr(2), mnarr(2), mxarr(2) & , avarr(3), mnarr(3), mxarr(3) & , avarr(4), mnarr(4), mxarr(4) & , avarr(5), mnarr(5), mxarr(5) & , avarr(6), mnarr(6), mxarr(6) & , avarr(7), mnarr(7), mxarr(7) call flush_and_sync(funit) call flush_and_sync(sunit) if (error_control) then write(eunit,efmt) step, time, dth, dte, maxval(errors(:)), errors(:) call flush_and_sync(eunit) end if end if #ifdef PROFILE call stop_timer(ims) #endif /* PROFILE */ !------------------------------------------------------------------------------- ! end subroutine store_integrals !=============================================================================== ! end module