messages.F90 Source File

MODULE messages


This file depends on

sourcefile~~messages.f90~~EfferentGraph sourcefile~messages.f90 messages.F90 sourcefile~constants_io.f90 constants_io.F90 sourcefile~messages.f90->sourcefile~constants_io.f90 sourcefile~precision.f90 precision.F90 sourcefile~messages.f90->sourcefile~precision.f90

Files dependent on this one

sourcefile~~messages.f90~~AfferentGraph sourcefile~messages.f90 messages.F90 sourcefile~apply_basis.f90 apply_basis.F90 sourcefile~apply_basis.f90->sourcefile~messages.f90 sourcefile~c_interop.f90 c_interop.F90 sourcefile~apply_basis.f90->sourcefile~c_interop.f90 sourcefile~tagarray_driver.f90 tagarray_driver.F90 sourcefile~apply_basis.f90->sourcefile~tagarray_driver.f90 sourcefile~types.f90 types.F90 sourcefile~apply_basis.f90->sourcefile~types.f90 sourcefile~basis_tools.f90 basis_tools.F90 sourcefile~basis_tools.f90->sourcefile~messages.f90 sourcefile~blas_wrap.f90 blas_wrap.F90 sourcefile~blas_wrap.f90->sourcefile~messages.f90 sourcefile~c_interop.f90->sourcefile~messages.f90 sourcefile~c_interop.f90->sourcefile~basis_tools.f90 sourcefile~c_interop.f90->sourcefile~tagarray_driver.f90 sourcefile~c_interop.f90->sourcefile~types.f90 sourcefile~dft.f90 dft.F90 sourcefile~dft.f90->sourcefile~messages.f90 sourcefile~dft.f90->sourcefile~basis_tools.f90 sourcefile~libxc.f90 libxc.F90 sourcefile~dft.f90->sourcefile~libxc.f90 sourcefile~mathlib.f90 mathlib.F90 sourcefile~dft.f90->sourcefile~mathlib.f90 sourcefile~dft.f90->sourcefile~tagarray_driver.f90 sourcefile~dft_fuzzycell.f90 dft_fuzzycell.F90 sourcefile~dft.f90->sourcefile~dft_fuzzycell.f90 sourcefile~dft_gridint_energy.f90 dft_gridint_energy.F90 sourcefile~dft.f90->sourcefile~dft_gridint_energy.f90 sourcefile~dft_gridint_grad.f90 dft_gridint_grad.F90 sourcefile~dft.f90->sourcefile~dft_gridint_grad.f90 sourcefile~dft_molgrid.f90 dft_molgrid.F90 sourcefile~dft.f90->sourcefile~dft_molgrid.f90 sourcefile~dft.f90->sourcefile~types.f90 sourcefile~eigen.f90 eigen.F90 sourcefile~eigen.f90->sourcefile~messages.f90 sourcefile~oqp_linalg.f90 oqp_linalg.F90 sourcefile~eigen.f90->sourcefile~oqp_linalg.f90 sourcefile~electric_moments.f90 electric_moments.F90 sourcefile~electric_moments.f90->sourcefile~messages.f90 sourcefile~electric_moments.f90->sourcefile~basis_tools.f90 sourcefile~electric_moments.f90->sourcefile~c_interop.f90 sourcefile~int1.f90 int1.F90 sourcefile~electric_moments.f90->sourcefile~int1.f90 sourcefile~electric_moments.f90->sourcefile~mathlib.f90 sourcefile~electric_moments.f90->sourcefile~tagarray_driver.f90 sourcefile~electric_moments.f90->sourcefile~types.f90 sourcefile~functionals.f90 functionals.F90 sourcefile~functionals.f90->sourcefile~messages.f90 sourcefile~get_basis_overlap.f90 get_basis_overlap.F90 sourcefile~get_basis_overlap.f90->sourcefile~messages.f90 sourcefile~get_basis_overlap.f90->sourcefile~basis_tools.f90 sourcefile~get_basis_overlap.f90->sourcefile~c_interop.f90 sourcefile~get_basis_overlap.f90->sourcefile~int1.f90 sourcefile~get_basis_overlap.f90->sourcefile~tagarray_driver.f90 sourcefile~get_basis_overlap.f90->sourcefile~types.f90 sourcefile~get_states_overlap.f90 get_states_overlap.F90 sourcefile~get_states_overlap.f90->sourcefile~messages.f90 sourcefile~get_states_overlap.f90->sourcefile~basis_tools.f90 sourcefile~get_states_overlap.f90->sourcefile~c_interop.f90 sourcefile~get_states_overlap.f90->sourcefile~tagarray_driver.f90 sourcefile~tdhf_mrsf_lib.f90 tdhf_mrsf_lib.F90 sourcefile~get_states_overlap.f90->sourcefile~tdhf_mrsf_lib.f90 sourcefile~get_states_overlap.f90->sourcefile~types.f90 sourcefile~grd1.f90 grd1.F90 sourcefile~grd1.f90->sourcefile~messages.f90 sourcefile~grd1.f90->sourcefile~basis_tools.f90 sourcefile~grd1.f90->sourcefile~mathlib.f90 sourcefile~grd1.f90->sourcefile~tagarray_driver.f90 sourcefile~ecp.f90 ecp.F90 sourcefile~grd1.f90->sourcefile~ecp.f90 sourcefile~mod_shell_tools.f90 mod_shell_tools.F90 sourcefile~grd1.f90->sourcefile~mod_shell_tools.f90 sourcefile~grd1.f90->sourcefile~types.f90 sourcefile~mod_1e_primitives.f90 mod_1e_primitives.F90 sourcefile~grd1.f90->sourcefile~mod_1e_primitives.f90 sourcefile~grd2.f90 grd2.F90 sourcefile~grd2.f90->sourcefile~messages.f90 sourcefile~grd2.f90->sourcefile~basis_tools.f90 sourcefile~int2.f90 int2.F90 sourcefile~grd2.f90->sourcefile~int2.f90 sourcefile~grd2.f90->sourcefile~tagarray_driver.f90 sourcefile~grd2.f90->sourcefile~dft_molgrid.f90 sourcefile~grd2_rys.f90 grd2_rys.F90 sourcefile~grd2.f90->sourcefile~grd2_rys.f90 sourcefile~int2_pairs.f90 int2_pairs.F90 sourcefile~grd2.f90->sourcefile~int2_pairs.f90 sourcefile~grd2.f90->sourcefile~types.f90 sourcefile~guess.f90 guess.F90 sourcefile~guess.f90->sourcefile~messages.f90 sourcefile~guess.f90->sourcefile~basis_tools.f90 sourcefile~guess.f90->sourcefile~eigen.f90 sourcefile~guess.f90->sourcefile~mathlib.f90 sourcefile~guess.f90->sourcefile~oqp_linalg.f90 sourcefile~guess.f90->sourcefile~types.f90 sourcefile~guess_hcore.f90 guess_hcore.F90 sourcefile~guess_hcore.f90->sourcefile~messages.f90 sourcefile~guess_hcore.f90->sourcefile~basis_tools.f90 sourcefile~guess_hcore.f90->sourcefile~c_interop.f90 sourcefile~guess_hcore.f90->sourcefile~guess.f90 sourcefile~guess_hcore.f90->sourcefile~mathlib.f90 sourcefile~printing.f90 printing.F90 sourcefile~guess_hcore.f90->sourcefile~printing.f90 sourcefile~guess_hcore.f90->sourcefile~tagarray_driver.f90 sourcefile~guess_hcore.f90->sourcefile~types.f90 sourcefile~guess_huckel.f90 guess_huckel.F90 sourcefile~guess_huckel.f90->sourcefile~messages.f90 sourcefile~guess_huckel.f90->sourcefile~basis_tools.f90 sourcefile~guess_huckel.f90->sourcefile~c_interop.f90 sourcefile~guess_huckel.f90->sourcefile~guess.f90 sourcefile~huckel.f90 huckel.F90 sourcefile~guess_huckel.f90->sourcefile~huckel.f90 sourcefile~guess_huckel.f90->sourcefile~mathlib.f90 sourcefile~guess_huckel.f90->sourcefile~printing.f90 sourcefile~guess_huckel.f90->sourcefile~tagarray_driver.f90 sourcefile~guess_huckel.f90->sourcefile~types.f90 sourcefile~guess_json.f90 guess_json.F90 sourcefile~guess_json.f90->sourcefile~messages.f90 sourcefile~guess_json.f90->sourcefile~basis_tools.f90 sourcefile~guess_json.f90->sourcefile~c_interop.f90 sourcefile~guess_json.f90->sourcefile~guess.f90 sourcefile~guess_json.f90->sourcefile~printing.f90 sourcefile~guess_json.f90->sourcefile~tagarray_driver.f90 sourcefile~guess_json.f90->sourcefile~types.f90 sourcefile~hf_energy.f90 hf_energy.f90 sourcefile~hf_energy.f90->sourcefile~messages.f90 sourcefile~hf_energy.f90->sourcefile~basis_tools.f90 sourcefile~hf_energy.f90->sourcefile~c_interop.f90 sourcefile~hf_energy.f90->sourcefile~dft.f90 sourcefile~hf_energy.f90->sourcefile~printing.f90 sourcefile~scf.f90 scf.F90 sourcefile~hf_energy.f90->sourcefile~scf.f90 sourcefile~hf_energy.f90->sourcefile~tagarray_driver.f90 sourcefile~hf_energy.f90->sourcefile~dft_molgrid.f90 sourcefile~hf_energy.f90->sourcefile~types.f90 sourcefile~hf_gradient.f90 hf_gradient.F90 sourcefile~hf_gradient.f90->sourcefile~messages.f90 sourcefile~hf_gradient.f90->sourcefile~basis_tools.f90 sourcefile~hf_gradient.f90->sourcefile~c_interop.f90 sourcefile~hf_gradient.f90->sourcefile~dft.f90 sourcefile~hf_gradient.f90->sourcefile~grd1.f90 sourcefile~hf_gradient.f90->sourcefile~grd2.f90 sourcefile~hf_gradient.f90->sourcefile~mathlib.f90 sourcefile~hf_gradient.f90->sourcefile~printing.f90 sourcefile~hf_gradient.f90->sourcefile~tagarray_driver.f90 sourcefile~hf_gradient.f90->sourcefile~dft_molgrid.f90 sourcefile~hf_gradient.f90->sourcefile~types.f90 sourcefile~huckel.f90->sourcefile~messages.f90 sourcefile~huckel.f90->sourcefile~basis_tools.f90 sourcefile~huckel.f90->sourcefile~eigen.f90 sourcefile~huckel.f90->sourcefile~guess.f90 sourcefile~huckel.f90->sourcefile~int1.f90 sourcefile~huckel.f90->sourcefile~mathlib.f90 sourcefile~huckel.f90->sourcefile~oqp_linalg.f90 sourcefile~huckel.f90->sourcefile~types.f90 sourcefile~int1.f90->sourcefile~messages.f90 sourcefile~int1.f90->sourcefile~basis_tools.f90 sourcefile~int1.f90->sourcefile~printing.f90 sourcefile~int1.f90->sourcefile~ecp.f90 sourcefile~int1.f90->sourcefile~mod_shell_tools.f90 sourcefile~int1.f90->sourcefile~mod_1e_primitives.f90 sourcefile~int1e.f90 int1e.F90 sourcefile~int1e.f90->sourcefile~messages.f90 sourcefile~int1e.f90->sourcefile~basis_tools.f90 sourcefile~int1e.f90->sourcefile~c_interop.f90 sourcefile~int1e.f90->sourcefile~int1.f90 sourcefile~int1e.f90->sourcefile~printing.f90 sourcefile~int1e.f90->sourcefile~tagarray_driver.f90 sourcefile~int1e.f90->sourcefile~types.f90 sourcefile~int2.f90->sourcefile~messages.f90 sourcefile~int2.f90->sourcefile~basis_tools.f90 sourcefile~int2.f90->sourcefile~tagarray_driver.f90 sourcefile~int2.f90->sourcefile~int2_pairs.f90 sourcefile~int_libint.f90 int_libint.F90 sourcefile~int2.f90->sourcefile~int_libint.f90 sourcefile~int_rotaxis.f90 int_rotaxis.F90 sourcefile~int2.f90->sourcefile~int_rotaxis.f90 sourcefile~int_rys.f90 int_rys.F90 sourcefile~int2.f90->sourcefile~int_rys.f90 sourcefile~int2.f90->sourcefile~types.f90 sourcefile~lapack_wrap.f90 lapack_wrap.F90 sourcefile~lapack_wrap.f90->sourcefile~messages.f90 sourcefile~lebedev.f90 lebedev.F90 sourcefile~lebedev.f90->sourcefile~messages.f90 sourcefile~libxc.f90->sourcefile~messages.f90 sourcefile~libxc.f90->sourcefile~functionals.f90 sourcefile~libxc.f90->sourcefile~types.f90 sourcefile~mathlib.f90->sourcefile~messages.f90 sourcefile~mathlib.f90->sourcefile~eigen.f90 sourcefile~mathlib.f90->sourcefile~oqp_linalg.f90 sourcefile~oqp_banner.f90 oqp_banner.F90 sourcefile~oqp_banner.f90->sourcefile~messages.f90 sourcefile~oqp_banner.f90->sourcefile~c_interop.f90 sourcefile~oqp_banner.f90->sourcefile~tagarray_driver.f90 sourcefile~oqp_banner.f90->sourcefile~types.f90 sourcefile~pcg.f90 pcg.F90 sourcefile~pcg.f90->sourcefile~messages.f90 sourcefile~population_analysis.f90 population_analysis.F90 sourcefile~population_analysis.f90->sourcefile~messages.f90 sourcefile~population_analysis.f90->sourcefile~basis_tools.f90 sourcefile~population_analysis.f90->sourcefile~c_interop.f90 sourcefile~population_analysis.f90->sourcefile~eigen.f90 sourcefile~population_analysis.f90->sourcefile~mathlib.f90 sourcefile~population_analysis.f90->sourcefile~tagarray_driver.f90 sourcefile~population_analysis.f90->sourcefile~oqp_linalg.f90 sourcefile~population_analysis.f90->sourcefile~types.f90 sourcefile~printing.f90->sourcefile~messages.f90 sourcefile~printing.f90->sourcefile~basis_tools.f90 sourcefile~printing.f90->sourcefile~tagarray_driver.f90 sourcefile~printing.f90->sourcefile~types.f90 sourcefile~resp.f90 resp.F90 sourcefile~resp.f90->sourcefile~messages.f90 sourcefile~resp.f90->sourcefile~basis_tools.f90 sourcefile~resp.f90->sourcefile~c_interop.f90 sourcefile~resp.f90->sourcefile~int1.f90 sourcefile~resp.f90->sourcefile~lebedev.f90 sourcefile~resp.f90->sourcefile~mathlib.f90 sourcefile~resp.f90->sourcefile~tagarray_driver.f90 sourcefile~resp.f90->sourcefile~oqp_linalg.f90 sourcefile~resp.f90->sourcefile~types.f90 sourcefile~scf.f90->sourcefile~messages.f90 sourcefile~scf.f90->sourcefile~basis_tools.f90 sourcefile~scf.f90->sourcefile~dft.f90 sourcefile~scf.f90->sourcefile~guess.f90 sourcefile~scf.f90->sourcefile~int2.f90 sourcefile~scf.f90->sourcefile~mathlib.f90 sourcefile~scf.f90->sourcefile~printing.f90 sourcefile~scf.f90->sourcefile~tagarray_driver.f90 sourcefile~scf.f90->sourcefile~dft_molgrid.f90 sourcefile~scf.f90->sourcefile~oqp_linalg.f90 sourcefile~scf_converger.f90 scf_converger.F90 sourcefile~scf.f90->sourcefile~scf_converger.f90 sourcefile~scf.f90->sourcefile~types.f90 sourcefile~tagarray_driver.f90->sourcefile~messages.f90 sourcefile~tdhf_energy.f90 tdhf_energy.F90 sourcefile~tdhf_energy.f90->sourcefile~messages.f90 sourcefile~tdhf_energy.f90->sourcefile~basis_tools.f90 sourcefile~tdhf_energy.f90->sourcefile~c_interop.f90 sourcefile~tdhf_energy.f90->sourcefile~dft.f90 sourcefile~tdhf_energy.f90->sourcefile~int1.f90 sourcefile~tdhf_energy.f90->sourcefile~int2.f90 sourcefile~tdhf_energy.f90->sourcefile~mathlib.f90 sourcefile~tdhf_energy.f90->sourcefile~printing.f90 sourcefile~tdhf_energy.f90->sourcefile~tagarray_driver.f90 sourcefile~dft_gridint_fxc.f90 dft_gridint_fxc.F90 sourcefile~tdhf_energy.f90->sourcefile~dft_gridint_fxc.f90 sourcefile~tdhf_energy.f90->sourcefile~dft_molgrid.f90 sourcefile~tdhf_lib.f90 tdhf_lib.F90 sourcefile~tdhf_energy.f90->sourcefile~tdhf_lib.f90 sourcefile~tdhf_energy.f90->sourcefile~types.f90 sourcefile~tdhf_gradient.f90 tdhf_gradient.F90 sourcefile~tdhf_gradient.f90->sourcefile~messages.f90 sourcefile~tdhf_gradient.f90->sourcefile~basis_tools.f90 sourcefile~tdhf_gradient.f90->sourcefile~c_interop.f90 sourcefile~tdhf_gradient.f90->sourcefile~dft.f90 sourcefile~tdhf_gradient.f90->sourcefile~grd1.f90 sourcefile~tdhf_gradient.f90->sourcefile~grd2.f90 sourcefile~tdhf_gradient.f90->sourcefile~mathlib.f90 sourcefile~tdhf_gradient.f90->sourcefile~printing.f90 sourcefile~tdhf_gradient.f90->sourcefile~tagarray_driver.f90 sourcefile~dft_gridint_tdxc_grad.f90 dft_gridint_tdxc_grad.F90 sourcefile~tdhf_gradient.f90->sourcefile~dft_gridint_tdxc_grad.f90 sourcefile~tdhf_gradient.f90->sourcefile~dft_molgrid.f90 sourcefile~tdhf_gradient.f90->sourcefile~tdhf_lib.f90 sourcefile~tdhf_gradient.f90->sourcefile~types.f90 sourcefile~tdhf_mrsf_energy.f90 tdhf_mrsf_energy.F90 sourcefile~tdhf_mrsf_energy.f90->sourcefile~messages.f90 sourcefile~tdhf_mrsf_energy.f90->sourcefile~basis_tools.f90 sourcefile~tdhf_mrsf_energy.f90->sourcefile~c_interop.f90 sourcefile~tdhf_mrsf_energy.f90->sourcefile~int2.f90 sourcefile~tdhf_mrsf_energy.f90->sourcefile~mathlib.f90 sourcefile~tdhf_mrsf_energy.f90->sourcefile~printing.f90 sourcefile~tdhf_mrsf_energy.f90->sourcefile~tagarray_driver.f90 sourcefile~tdhf_mrsf_energy.f90->sourcefile~tdhf_mrsf_lib.f90 sourcefile~tdhf_sf_lib.f90 tdhf_sf_lib.F90 sourcefile~tdhf_mrsf_energy.f90->sourcefile~tdhf_sf_lib.f90 sourcefile~tdhf_mrsf_energy.f90->sourcefile~oqp_linalg.f90 sourcefile~tdhf_mrsf_energy.f90->sourcefile~tdhf_lib.f90 sourcefile~tdhf_mrsf_energy.f90->sourcefile~types.f90 sourcefile~tdhf_mrsf_lib.f90->sourcefile~messages.f90 sourcefile~tdhf_mrsf_lib.f90->sourcefile~basis_tools.f90 sourcefile~tdhf_mrsf_lib.f90->sourcefile~int2.f90 sourcefile~tdhf_mrsf_lib.f90->sourcefile~oqp_linalg.f90 sourcefile~tdhf_mrsf_lib.f90->sourcefile~types.f90 sourcefile~tdhf_sf_energy.f90 tdhf_sf_energy.F90 sourcefile~tdhf_sf_energy.f90->sourcefile~messages.f90 sourcefile~tdhf_sf_energy.f90->sourcefile~basis_tools.f90 sourcefile~tdhf_sf_energy.f90->sourcefile~c_interop.f90 sourcefile~tdhf_sf_energy.f90->sourcefile~int2.f90 sourcefile~tdhf_sf_energy.f90->sourcefile~mathlib.f90 sourcefile~tdhf_sf_energy.f90->sourcefile~printing.f90 sourcefile~tdhf_sf_energy.f90->sourcefile~tagarray_driver.f90 sourcefile~tdhf_sf_energy.f90->sourcefile~tdhf_sf_lib.f90 sourcefile~tdhf_sf_energy.f90->sourcefile~oqp_linalg.f90 sourcefile~tdhf_sf_energy.f90->sourcefile~tdhf_lib.f90 sourcefile~tdhf_sf_energy.f90->sourcefile~types.f90 sourcefile~tdhf_sf_lib.f90->sourcefile~messages.f90 sourcefile~tdhf_sf_lib.f90->sourcefile~basis_tools.f90 sourcefile~tdhf_sf_lib.f90->sourcefile~int1.f90 sourcefile~tdhf_sf_lib.f90->sourcefile~mathlib.f90 sourcefile~tdhf_sf_lib.f90->sourcefile~oqp_linalg.f90 sourcefile~tdhf_sf_lib.f90->sourcefile~tdhf_lib.f90 sourcefile~tdhf_sf_lib.f90->sourcefile~types.f90 sourcefile~tdhf_z_vector.f90 tdhf_z_vector.F90 sourcefile~tdhf_z_vector.f90->sourcefile~messages.f90 sourcefile~tdhf_z_vector.f90->sourcefile~basis_tools.f90 sourcefile~tdhf_z_vector.f90->sourcefile~c_interop.f90 sourcefile~tdhf_z_vector.f90->sourcefile~dft.f90 sourcefile~tdhf_z_vector.f90->sourcefile~int2.f90 sourcefile~tdhf_z_vector.f90->sourcefile~mathlib.f90 sourcefile~tdhf_z_vector.f90->sourcefile~pcg.f90 sourcefile~tdhf_z_vector.f90->sourcefile~printing.f90 sourcefile~tdhf_z_vector.f90->sourcefile~tagarray_driver.f90 sourcefile~tdhf_z_vector.f90->sourcefile~tdhf_sf_lib.f90 sourcefile~tdhf_z_vector.f90->sourcefile~dft_gridint_fxc.f90 sourcefile~dft_gridint_gxc.f90 dft_gridint_gxc.F90 sourcefile~tdhf_z_vector.f90->sourcefile~dft_gridint_gxc.f90 sourcefile~tdhf_z_vector.f90->sourcefile~dft_molgrid.f90 sourcefile~tdhf_z_vector.f90->sourcefile~oqp_linalg.f90 sourcefile~tdhf_z_vector.f90->sourcefile~tdhf_lib.f90 sourcefile~tdhf_z_vector.f90->sourcefile~types.f90 sourcefile~dft_fuzzycell.f90->sourcefile~basis_tools.f90 sourcefile~dft_fuzzycell.f90->sourcefile~dft_molgrid.f90 sourcefile~dft_gridint.f90 dft_gridint.F90 sourcefile~dft_gridint.f90->sourcefile~basis_tools.f90 sourcefile~dft_gridint.f90->sourcefile~functionals.f90 sourcefile~dft_gridint.f90->sourcefile~dft_molgrid.f90 sourcefile~dft_xc_libxc.f90 dft_xc_libxc.F90 sourcefile~dft_gridint.f90->sourcefile~dft_xc_libxc.f90 sourcefile~dft_gridint.f90->sourcefile~oqp_linalg.f90 sourcefile~dft_gridint_energy.f90->sourcefile~basis_tools.f90 sourcefile~dft_gridint_energy.f90->sourcefile~dft_gridint.f90 sourcefile~dft_gridint_energy.f90->sourcefile~dft_molgrid.f90 sourcefile~dft_gridint_energy.f90->sourcefile~oqp_linalg.f90 sourcefile~dft_gridint_energy.f90->sourcefile~types.f90 sourcefile~dft_gridint_fxc.f90->sourcefile~basis_tools.f90 sourcefile~dft_gridint_fxc.f90->sourcefile~mathlib.f90 sourcefile~dft_gridint_fxc.f90->sourcefile~dft_gridint.f90 sourcefile~dft_gridint_fxc.f90->sourcefile~dft_molgrid.f90 sourcefile~dft_gridint_fxc.f90->sourcefile~oqp_linalg.f90 sourcefile~dft_gridint_fxc.f90->sourcefile~types.f90 sourcefile~dft_gridint_grad.f90->sourcefile~basis_tools.f90 sourcefile~dft_gridint_grad.f90->sourcefile~dft_gridint.f90 sourcefile~dft_gridint_grad.f90->sourcefile~dft_molgrid.f90 sourcefile~dft_gridint_grad.f90->sourcefile~types.f90 sourcefile~dft_gridint_gxc.f90->sourcefile~basis_tools.f90 sourcefile~dft_gridint_gxc.f90->sourcefile~mathlib.f90 sourcefile~dft_gridint_gxc.f90->sourcefile~dft_gridint.f90 sourcefile~dft_gridint_gxc.f90->sourcefile~dft_gridint_fxc.f90 sourcefile~dft_gridint_gxc.f90->sourcefile~dft_molgrid.f90 sourcefile~dft_gridint_gxc.f90->sourcefile~oqp_linalg.f90 sourcefile~dft_gridint_gxc.f90->sourcefile~types.f90 sourcefile~dft_gridint_tdxc_grad.f90->sourcefile~basis_tools.f90 sourcefile~dft_gridint_tdxc_grad.f90->sourcefile~dft_gridint.f90 sourcefile~dft_gridint_tdxc_grad.f90->sourcefile~dft_molgrid.f90 sourcefile~dft_gridint_tdxc_grad.f90->sourcefile~types.f90 sourcefile~dft_molgrid.f90->sourcefile~lebedev.f90 sourcefile~dft_xc_libxc.f90->sourcefile~functionals.f90 sourcefile~dft_xclib.f90 dft_xclib.F90 sourcefile~dft_xc_libxc.f90->sourcefile~dft_xclib.f90 sourcefile~dft_xclib.f90->sourcefile~functionals.f90 sourcefile~ecp.f90->sourcefile~basis_tools.f90 sourcefile~grd2_rys.f90->sourcefile~basis_tools.f90 sourcefile~grd2_rys.f90->sourcefile~int2_pairs.f90 sourcefile~int2_pairs.f90->sourcefile~basis_tools.f90 sourcefile~int_libint.f90->sourcefile~basis_tools.f90 sourcefile~int_libint.f90->sourcefile~int2_pairs.f90 sourcefile~int_rotaxis.f90->sourcefile~basis_tools.f90 sourcefile~int_rotaxis.f90->sourcefile~int2_pairs.f90 sourcefile~int_rys.f90->sourcefile~basis_tools.f90 sourcefile~int_rys.f90->sourcefile~int2_pairs.f90 sourcefile~mod_shell_tools.f90->sourcefile~basis_tools.f90 sourcefile~oqp_linalg.f90->sourcefile~blas_wrap.f90 sourcefile~oqp_linalg.f90->sourcefile~lapack_wrap.f90 sourcefile~scf_converger.f90->sourcefile~mathlib.f90 sourcefile~scf_converger.f90->sourcefile~printing.f90 sourcefile~scf_converger.f90->sourcefile~oqp_linalg.f90 sourcefile~tdhf_lib.f90->sourcefile~basis_tools.f90 sourcefile~tdhf_lib.f90->sourcefile~eigen.f90 sourcefile~tdhf_lib.f90->sourcefile~int2.f90 sourcefile~tdhf_lib.f90->sourcefile~mathlib.f90 sourcefile~tdhf_lib.f90->sourcefile~oqp_linalg.f90 sourcefile~tdhf_lib.f90->sourcefile~types.f90 sourcefile~types.f90->sourcefile~basis_tools.f90 sourcefile~types.f90->sourcefile~functionals.f90 sourcefile~mod_1e_primitives.f90->sourcefile~mod_shell_tools.f90

Source Code

!*MODULE messages
!> @brief   This module provides routines where the output is
!> @details Mostly, this file is needed for simplifying of usage
!>            the LibXC interface in different software
!>          For GAMESS(US), this file can be expanded for other messages
!>            For example, aborting with printing custom message
!> @author  Igor S. Gerasimov
!> @date    July, 2021 - Initial release -
!> @todo    Remove it
!> @params  WITH_ABORT    - logical key for stopmode; stop will be
!> @params  WITHOUT_ABORT - logical key for stopmode; stop will not be
module messages
  use precision, only: dp
#ifdef OQP
  use io_constants, only: write_unit => iw
#else
  use comm_IOFILE, only: write_unit => IW
  use comm_PAR, only: master_worker => MASWRK
#endif
  implicit none
  private
  logical, parameter :: WITH_ABORT = .true.
  logical, parameter :: WITHOUT_ABORT = .false.
  interface show_message
    module procedure show_message_text, &
      show_message_with_integer, &
      show_message_with_double, &
      show_message_with_double_and_text, &
      show_message_with_integer_and_text, &
      show_message_with_keys
  end interface show_message
  public show_message, with_abort, without_abort
contains
  !> @brief   Print simple message
  !> @author  Igor S. Gerasimov
  !> @date    July, 2021 - Initial release -
  !> @params  message  (in)           - message for displaying
  !> @params  stopmode (in, optional) - is aborting required?
  subroutine show_message_text(message, stopmode)
    character(len=*), intent(in) :: message
    logical, intent(in), optional :: stopmode
    logical :: stopmode_
    if (.not. present(stopmode)) then
      stopmode_ = WITHOUT_ABORT
    else
      stopmode_ = stopmode
    end if
#ifndef OQP
    if (master_worker) &
#endif
      write (write_unit, "(A)") message
    call abort(stopmode_)
  end subroutine show_message_text
  !> @brief   Print simple message
  !> @details write( ,format) message, value
  !> @author  Igor S. Gerasimov
  !> @date    July, 2021 - Initial release -
  !> @params  format   (in)           - format for displaying
  !> @params  message  (in)           - message for displaying
  !> @params  value    (in)           - value for displaying
  !> @params  stopmode (in, optional) - is aborting required?
  subroutine show_message_with_integer(format, message, value, stopmode)
    character(len=*), intent(in) :: format
    character(len=*), intent(in) :: message
    integer, intent(in) :: value
    logical, intent(in), optional :: stopmode
    logical :: stopmode_
    if (.not. present(stopmode)) then
      stopmode_ = WITHOUT_ABORT
    else
      stopmode_ = stopmode
    end if
#ifndef OQP
    if (master_worker) &
#endif
      write (write_unit, format) message, value
    call abort(stopmode_)
  end subroutine show_message_with_integer
  !> @brief   Print simple message
  !> @details write( ,format) message, value
  !> @author  Igor S. Gerasimov
  !> @date    July, 2021 - Initial release -
  !> @params  format   (in)           - format for displaying
  !> @params  message  (in)           - message for displaying
  !> @params  value    (in)           - value for displaying
  !> @params  stopmode (in, optional) - is aborting required?
  subroutine show_message_with_double(format, message, value, stopmode)
    character(len=*), intent(in) :: format
    character(len=*), intent(in) :: message
    real(dp), intent(in) :: value
    logical, intent(in), optional :: stopmode
    logical :: stopmode_
    if (.not. present(stopmode)) then
      stopmode_ = WITHOUT_ABORT
    else
      stopmode_ = stopmode
    end if
#ifndef OQP
    if (master_worker) &
#endif
      write (write_unit, format) message, value
    call abort(stopmode_)
  end subroutine show_message_with_double
  !> @brief   Print simple message
  !> @details write( ,format) message1, value, message2
  !> @author  Igor S. Gerasimov
  !> @date    July, 2021 - Initial release -
  !> @params  format   (in)           - format for displaying
  !> @params  message1 (in)           - message for displaying
  !> @params  value    (in)           - value for displaying
  !> @params  message2 (in)           - message for displaying
  !> @params  stopmode (in, optional) - is aborting required?
  subroutine show_message_with_double_and_text(format, message1, value, message2, stopmode)
    character(len=*), intent(in) :: format
    character(len=*), intent(in) :: message1, message2
    real(kind=dp), intent(in) :: value
    logical, intent(in), optional :: stopmode
    logical :: stopmode_
    if (.not. present(stopmode)) then
      stopmode_ = WITHOUT_ABORT
    else
      stopmode_ = stopmode
    end if
#ifndef OQP
    if (master_worker) &
#endif
      write (write_unit, format) message1, value, message2
    call abort(stopmode_)
  end subroutine show_message_with_double_and_text
  !> @brief   Print simple message
  !> @details write( ,format) message1, value, message2
  !> @author  Igor S. Gerasimov
  !> @date    July, 2021 - Initial release -
  !> @params  format   (in)           - format for displaying
  !> @params  message1 (in)           - message for displaying
  !> @params  value    (in)           - value for displaying
  !> @params  message2 (in)           - message for displaying
  !> @params  stopmode (in, optional) - is aborting required?
  subroutine show_message_with_integer_and_text(format, message1, value, message2, stopmode)
    character(len=*), intent(in) :: format
    character(len=*), intent(in) :: message1, message2
    integer, intent(in) :: value
    logical, intent(in), optional :: stopmode
    logical :: stopmode_
    if (.not. present(stopmode)) then
      stopmode_ = WITHOUT_ABORT
    else
      stopmode_ = stopmode
    end if
#ifndef OQP
    if (master_worker) &
#endif
      write (write_unit, format) message1, value, message2
    call abort(stopmode_)
  end subroutine show_message_with_integer_and_text
  !> @brief   Print simple message
  !> @details write( ,'(A)') message
  !>          write( ,format) keys
  !> @author  Igor S. Gerasimov
  !> @date    July, 2021 - Initial release -
  !> @params  message  (in)           - message for displaying
  !> @params  format   (in)           - format for displaying
  !> @params  keys     (in)           - keys for displaying
  !> @params  stopmode (in, optional) - is aborting required?
  subroutine show_message_with_keys(message, format, keys, stopmode)
    character(len=*), intent(in) :: message
    character(len=*), intent(in) :: format
    character(len=*), intent(in), optional :: keys(:)
    logical, intent(in), optional :: stopmode
    logical :: stopmode_
    if (.not. present(stopmode)) then
      stopmode_ = WITHOUT_ABORT
    else
      stopmode_ = stopmode
    end if
#ifndef OQP
    if (master_worker) then
#endif
      write (write_unit, "(A)") message
      write (write_unit, format) keys
#ifndef OQP
    end if
#endif
    call abort(stopmode_)
  end subroutine show_message_with_keys
  subroutine abort(stopmode)
    logical, intent(in) :: stopmode
    flush(write_unit)
    if (stopmode) then
#ifdef OQP
#ifdef __GFORTRAN__
      call backtrace
#endif
      error stop "See above"
#else
      call abrt
#endif
    end if
  end subroutine abort
end module messages