27  use, 
intrinsic :: iso_fortran_env
 
   91  integer, 
parameter :: PCM_DIM_SPACE = 3
 
   95    logical, 
public                  :: run_pcm
 
   96    integer, 
public                  :: tdlevel
 
   98    integer, 
public                  :: n_tesserae
 
   99    type(pcm_sphere_t),  
allocatable :: spheres(:)
 
  100    type(pcm_tessera_t), 
allocatable, 
public :: tess(:)
 
  101    real(real64)                     :: scale_r
 
  102    real(real64), 
allocatable               :: matrix(:,:)
 
  103    real(real64), 
allocatable               :: matrix_d(:,:)
 
  104    real(real64), 
allocatable               :: matrix_lf(:,:)
 
  105    real(real64), 
allocatable               :: matrix_lf_d(:,:)
 
  106    real(real64), 
allocatable, 
public       :: q_e(:)
 
  107    real(real64), 
allocatable               :: q_n(:)
 
  108    real(real64), 
allocatable, 
public       :: q_e_in(:)
 
  109    real(real64), 
allocatable               :: rho_e(:)
 
  110    real(real64), 
allocatable               :: rho_n(:)
 
  111    real(real64)                     :: qtot_e
 
  112    real(real64)                     :: qtot_n
 
  113    real(real64)                     :: qtot_e_in
 
  114    real(real64)                     :: q_e_nominal
 
  115    real(real64)                     :: q_n_nominal
 
  116    logical                          :: renorm_charges
 
  117    real(real64)                     :: q_tot_tol
 
  118    real(real64)                     :: deltaQ_e
 
  119    real(real64)                     :: deltaQ_n
 
  120    real(real64), 
allocatable               :: v_e(:)
 
  121    real(real64), 
allocatable               :: v_n(:)
 
  122    real(real64), 
allocatable, 
public       :: v_e_rs(:)
 
  123    real(real64), 
allocatable, 
public       :: v_n_rs(:)
 
  124    real(real64), 
allocatable               :: q_ext(:)
 
  125    real(real64), 
allocatable               :: q_ext_in(:)
 
  126    real(real64), 
allocatable               :: rho_ext(:)
 
  127    real(real64)                     :: qtot_ext
 
  128    real(real64)                     :: qtot_ext_in
 
  129    real(real64), 
allocatable               :: v_ext(:)
 
  130    real(real64), 
allocatable, 
public       :: v_ext_rs(:)
 
  131    real(real64), 
allocatable               :: q_kick(:)
 
  132    real(real64), 
allocatable               :: rho_kick(:)
 
  133    real(real64)                     :: qtot_kick
 
  134    real(real64), 
allocatable               :: v_kick(:)
 
  135    real(real64), 
allocatable, 
public       :: v_kick_rs(:)
 
  136    real(real64), 
public                    :: epsilon_0
 
  137    real(real64), 
public                    :: epsilon_infty
 
  138    integer, 
public                  :: which_eps
 
  139    type(debye_param_t)              :: deb
 
  140    type(drude_param_t)              :: drl
 
  141    logical, 
public                  :: localf
 
  142    logical, 
public                  :: solute
 
  143    logical                          :: kick_is_present
 
  145    logical, 
public                  :: kick_like
 
  146    integer                          :: initial_asc
 
  147    real(real64)                     :: gaussian_width
 
  149    integer, 
public                  :: counter
 
  150    character(len=80)                :: input_cavity
 
  151    integer                          :: update_iter
 
  152    integer, 
public                  :: iter
 
  153    integer                          :: calc_method
 
  155    real(real64), 
public                    :: dt
 
  158    type(namespace_t), 
pointer       :: namespace
 
  159    type(space_t)                    :: space
 
  168    type(debye_param_t)  :: deb
 
  169    type(drude_param_t)  :: drl
 
  172  real(real64), 
allocatable :: s_mat_act(:,:)
 
  173  real(real64), 
allocatable :: d_mat_act(:,:)
 
  174  real(real64), 
allocatable :: Sigma(:,:)
 
  175  real(real64), 
allocatable :: Delta(:,:)
 
  180  integer, 
parameter :: &
 
  185  integer, 
parameter, 
public :: &
 
  186    PCM_CALC_DIRECT  = 1,       &
 
  189  integer, 
parameter ::    &
 
  200  subroutine pcm_init(pcm, namespace, space, ions, grid, qtot, val_charge, external_potentials_present, kick_present)
 
  206    real(real64),              
intent(in)  :: qtot
 
  207    real(real64),              
intent(in)  :: val_charge
 
  208    logical,                   
intent(in)  :: external_potentials_present
 
  209    logical,                   
intent(in)  :: kick_present
 
  211    integer :: ia, ii, itess, jtess, pcm_vdw_type, subdivider
 
  212    integer :: cav_unit_test, iunit, pcmmat_unit
 
  213    integer :: pcmmat_gamess_unit, cav_gamess_unit
 
  214    real(real64)   :: min_distance
 
  216    integer, 
parameter :: mxts = 10000
 
  218    real(real64) :: default_value
 
  219    real(real64) :: vdw_radius
 
  224    logical :: add_spheres_h
 
  225    logical :: changed_default_nn
 
  227    integer :: default_nn
 
  228    real(real64)   :: max_area
 
  232    pcm%kick_is_present = kick_present
 
  237    pcm%kick_like = .false.
 
  239    pcm%namespace => namespace
 
  254    call parse_variable(namespace, 
'PCMCalculation', .false., pcm%run_pcm)
 
  255    if (pcm%run_pcm) 
then 
  257      if (pcm%space%dim /= pcm_dim_space) 
then 
  258        message(1) = 
"PCM is only available for 3d calculations" 
  261      select type (box => grid%box)
 
  264        message(1) = 
"PCM is only available for BoxShape = minimum" 
  288    select case (pcm_vdw_type)
 
  290      default_value = 1.2_real64
 
  292      default_value = 
m_one 
  303    call parse_variable(namespace, 
'PCMRadiusScaling', default_value, pcm%scale_r)
 
  327    if (pcm%tdlevel /= 
pcm_td_eq .and. (.not. pcm%run_pcm)) 
then 
  328      call messages_write(
'Sorry, you have set PCMTDLevel /= eq, but PCMCalculation = no.')
 
  330      call messages_write(
'To spare you some time, Octopus will proceed as if PCMCalculation = yes.')
 
  352    call parse_variable(namespace, 
'PCMDynamicEpsilon', pcm%epsilon_0, pcm%epsilon_infty)
 
  375      call messages_write(
'Sorry, only Debye or Drude-Lorentz dielectric models are available.')
 
  377      call messages_write(
'To spare you some time, Octopus will proceed with the default choice (Debye).')
 
  379      call messages_write(
'You may change PCMEpsilonModel value for a Drude-Lorentz run.')
 
  384      call messages_write(
'Sorry, inertial/dynamic polarization splitting scheme for TD-PCM')
 
  387      call messages_write(
'require both static and dynamic dielectric constants,')
 
  390      call messages_write(
'Octopus will run using TD-PCM version in equilibrium')
 
  410    call parse_variable(namespace, 
'PCMEoMInitialCharges', 0, pcm%initial_asc)
 
  413    if (pcm%initial_asc /= 0) 
then 
  417        call messages_write(
'Sorry, initial polarization charges can only be read')
 
  421        call messages_write(
'Octopus will proceed as if PCMEoMInitialCharges = 0.')
 
  428    pcm%deb%eps_0 = pcm%epsilon_0
 
  429    pcm%deb%eps_d = pcm%epsilon_infty
 
  446      (abs(pcm%deb%tau) <= 
m_epsilon .or. 
is_close(pcm%deb%eps_0, pcm%deb%eps_d))) 
then 
  448      call messages_write(
'but you have not included all required Debye model parameters.')
 
  450      call messages_write(
'You need PCMEpsilonStatic, PCMEpsilonDynamic')
 
  451      call messages_write(
'and PCMDebyeRelaxTime for an EOM TD-PCM run.')
 
  453      call messages_write(
'Octopus will run using TD-PCM version in equilibrium')
 
  461        message(1) = 
"PCMEpsilonStatic = 1 is incompatible with a Drude-Lorentz EOM-PCM run." 
  480      call messages_write(
'but this is incompatible with a Drude-Lorentz EOM-PCM run.')
 
  483        call messages_write(
'Octopus will run using the default value of PCMDrudeLOmega.')
 
  487        message(1) = 
"PCMEpsilonStatic = 1 is incompatible with a Drude-Lorentz EOM-PCM run." 
  504    pcm%drl%aa = (pcm%epsilon_0 - 
m_one)*pcm%drl%w0**2
 
  516    call parse_variable(namespace, 
'PCMLocalField', .false., pcm%localf)
 
  519    if (pcm%localf .and. ((.not. external_potentials_present) .and. (.not. pcm%kick_is_present))) 
then 
  520      message(1) = 
"Sorry, you have set PCMLocalField = yes, but you have not included any external potentials." 
  536    if (pcm%run_pcm .and. (.not. pcm%solute)) 
then 
  537      call messages_write(
'N.B. This PCM run do not consider the polarization effects due to the solute.')
 
  539      if (.not. pcm%localf) 
then 
  540        message(1) = 
"You have activated a PCM run without polarization effects. Octopus will halt." 
  559    if (pcm%kick_like .and. (.not. pcm%run_pcm)) 
then 
  560      message(1) = 
"PCMKick option can only be activated when PCMCalculation = yes. Octopus will halt." 
  564    if (pcm%kick_like .and. (.not. pcm%localf)) 
then 
  565      message(1) = 
"PCMKick option can only be activated when a PCMLocalField = yes. Octopus will halt." 
  569    if (pcm%kick_like .and. (.not. pcm%kick_is_present)) 
then 
  570      message(1) = 
"Sorry, you have set PCMKick = yes, but you have not included any kick." 
  574    if (pcm%kick_is_present .and. pcm%run_pcm .and. (.not. pcm%localf)) 
then 
  575      message(1) = 
"You have set up a PCM calculation with a kick without local field effects." 
  576      message(2) = 
"Please, reconsider if you want the kick to be relevant for the PCM run." 
  587    call parse_variable(namespace, 
'PCMUpdateIter', 1, pcm%update_iter)
 
  611    call parse_variable(namespace, 
'PCMRenormCharges', .false., pcm%renorm_charges)
 
  627    if (pcm%renorm_charges) 
then 
  628      message(1) = 
"Info: Polarization charges will be renormalized" 
  629      message(2) = 
"      if |Q_tot_PCM - Q_M| > PCMQtotTol" 
  645    if (abs(pcm%gaussian_width) <= 
m_epsilon) 
then 
  646      message(1) = 
"Info: PCM potential will be defined in terms of polarization point charges" 
  649      message(1) = 
"Info: PCM potential is regularized to avoid Coulomb singularity" 
  671    if (pcm%input_cavity == 
'') 
then 
  680      call parse_variable(namespace, 
'PCMSpheresOnH', .false., add_spheres_h)
 
  684      do ia = 1, ions%natoms
 
  685        if ((.not. add_spheres_h) .and. ions%atom(ia)%label == 
'H') cycle
 
  686        pcm%n_spheres = pcm%n_spheres + 1 
 
  689      safe_allocate(pcm%spheres(1:pcm%n_spheres))
 
  696      do ia = 1, ions%natoms
 
  697        if ((.not. add_spheres_h) .and. ions%atom(ia)%label == 
'H') cycle
 
  698        pcm%n_spheres = pcm%n_spheres + 1
 
  701        pcm%spheres(pcm%n_spheres)%x = ions%pos(1, ia)
 
  702        pcm%spheres(pcm%n_spheres)%y = ions%pos(2, ia)
 
  703        pcm%spheres(pcm%n_spheres)%z = ions%pos(3, ia)
 
  706        pcm%spheres(pcm%n_spheres)%r = vdw_radius*pcm%scale_r
 
  710        pcm%info_unit = 
io_open(
pcm_dir//
'pcm_info.out', pcm%namespace, action=
'write')
 
  712        write(pcm%info_unit, 
'(A35)') 
'# Configuration: Molecule + Solvent' 
  713        write(pcm%info_unit, 
'(A35)') 
'# ---------------------------------' 
  714        write(pcm%info_unit, 
'(A21,F12.3)') 
'# Epsilon(Solvent) = ', pcm%epsilon_0
 
  715        write(pcm%info_unit, 
'(A1)')
'#' 
  716        write(pcm%info_unit, 
'(A35,I4)') 
'# Number of interlocking spheres = ', pcm%n_spheres
 
  717        write(pcm%info_unit, 
'(A1)')
'#' 
  719        write(pcm%info_unit, 
'(A8,3X,A7,8X,A26,20X,A10)') 
'# SPHERE', 
'ELEMENT', 
'CENTER  (X,Y,Z) (A)', 
'RADIUS (A)' 
  720        write(pcm%info_unit, 
'(A8,3X,A7,4X,A43,7X,A10)') 
'# ------', 
'-------', &
 
  721          '-------------------------------------------', 
'----------' 
  725      do ia = 1, ions%natoms
 
  726        if ((.not. add_spheres_h) .and. ions%atom(ia)%label == 
'H') cycle
 
  727        pcm%n_spheres = pcm%n_spheres + 1
 
  729          write(pcm%info_unit,
'(A1,2X,I3,7X,A2,3X,F14.8,2X,F14.8,2X,F14.8,4X,F14.8)')
'#', pcm%n_spheres, &
 
  730            ions%atom(ia)%label,          &
 
  731            ions%pos(1, ia)*
p_a_b,     &
 
  732            ions%pos(2, ia)*
p_a_b,     &
 
  733            ions%pos(3, ia)*
p_a_b,     &
 
  734            pcm%spheres(pcm%n_spheres)%r*
p_a_b 
  746      call parse_variable(namespace, 
'PCMTessSubdivider', 1, subdivider)
 
  748      safe_allocate(dum2(1:subdivider*
n_tess_sphere*pcm%n_spheres))
 
  758      call parse_variable(namespace, 
'PCMTessMinDistance', 0.1_real64, min_distance)
 
  762      call cav_gen(subdivider, min_distance, pcm%n_spheres, pcm%spheres, pcm%n_tesserae, dum2, pcm%info_unit)
 
  764      safe_allocate(pcm%tess(1:pcm%n_tesserae))
 
  766      do ia=1, pcm%n_tesserae
 
  767        pcm%tess(ia)%point    = 
m_zero 
  768        pcm%tess(ia)%area     = 
m_zero 
  769        pcm%tess(ia)%r_sphere = 
m_zero 
  770        pcm%tess(ia)%normal   = 
m_zero 
  773      pcm%tess = dum2(1:pcm%n_tesserae)
 
  775      safe_deallocate_a(dum2)
 
  777      message(1) = 
"Info: van der Waals surface has been calculated" 
  783      iunit = 
io_open(trim(pcm%input_cavity), pcm%namespace, action=
'read', status=
'old')
 
  784      read(iunit,*) pcm%n_tesserae
 
  786      if (pcm%n_tesserae > mxts) 
then 
  787        write(
message(1),
'(a,I5,a,I5)') 
"total number of tesserae", pcm%n_tesserae, 
">", mxts
 
  791      safe_allocate(pcm%tess(1:pcm%n_tesserae))
 
  793      do ia = 1, pcm%n_tesserae
 
  794        pcm%tess(ia)%point    = 
m_zero 
  795        pcm%tess(ia)%area     = 
m_zero 
  796        pcm%tess(ia)%r_sphere = 
m_zero 
  797        pcm%tess(ia)%normal   = 
m_zero 
  800      do ia = 1, pcm%n_tesserae
 
  801        read(iunit,*) pcm%tess(ia)%point(1)
 
  804      do ia = 1, pcm%n_tesserae
 
  805        read(iunit,*) pcm%tess(ia)%point(2)
 
  808      do ia = 1, pcm%n_tesserae
 
  809        read(iunit,*) pcm%tess(ia)%point(3)
 
  812      do ia = 1, pcm%n_tesserae
 
  813        read(iunit,*) pcm%tess(ia)%area
 
  816      do ia = 1, pcm%n_tesserae
 
  817        read(iunit,*) pcm%tess(ia)%r_sphere
 
  820      do ia = 1, pcm%n_tesserae
 
  821        read(iunit,*) pcm%tess(ia)%normal
 
  825      message(1) = 
"Info: van der Waals surface has been read from " // trim(pcm%input_cavity)
 
  830      cav_unit_test = 
io_open(
pcm_dir//
'cavity_mol.xyz', pcm%namespace, action=
'write')
 
  832      write (cav_unit_test,
'(2X,I4)') pcm%n_tesserae + ions%natoms
 
  833      write (cav_unit_test,
'(2X)')
 
  835      do ia = 1, pcm%n_tesserae
 
  836        write(cav_unit_test,
'(2X,A2,3X,4f15.8,3X,4f15.8,3X,4f15.8)') 
'H', pcm%tess(ia)%point*
p_a_b 
  839      do ia = 1, ions%natoms
 
  840        write(cav_unit_test,
'(2X,A2,3X,4f15.8,3X,4f15.8,3X,4f15.8)') ions%atom(ia)%label,      &
 
  841          ions%pos(:, ia)*
p_a_b 
  846      write(pcm%info_unit,
'(A1)')
'#' 
  848        write(pcm%info_unit,
'(A1,4X,A4,14X,A4,21X,A4,21X,A4,21X,A4,21X,A7,18X,A7,18X,A8,17X,A5,20X,A8,17X,A5,20X,A8,17X,A5)') &
 
  849          '#',
'iter', 
'E_ee', 
'E_en', 
'E_nn', 
'E_ne', 
'E_e_ext', 
'E_n_ext', 
'E_M-solv', &
 
  850          'Q_pol^e', 
'deltaQ^e', 
'Q_pol^n', 
'deltaQ^n', 
'Q_pol^ext' 
  852        write(pcm%info_unit,
'(A1,4X,A4,14X,A4,21X,A4,21X,A4,21X,A4,21X,A8,17X,A5,20X,A8,17X,A5,20X, A8)') &
 
  853          '#',
'iter', 
'E_ee', 
'E_en', 
'E_nn', 
'E_ne', 
'E_M-solv', 
'Q_pol^e', 
'deltaQ^e', 
'Q_pol^n', 
'deltaQ^n' 
  860      cav_gamess_unit = 
io_open(
pcm_dir//
'geom_cavity_gamess.out', pcm%namespace, action=
'write')
 
  862      write(cav_gamess_unit,*) pcm%n_tesserae
 
  864      do ia = 1, pcm%n_tesserae
 
  865        write(cav_gamess_unit,*) pcm%tess(ia)%point(1)
 
  868      do ia = 1, pcm%n_tesserae
 
  869        write(cav_gamess_unit,*) pcm%tess(ia)%point(2)
 
  872      do ia = 1, pcm%n_tesserae
 
  873        write(cav_gamess_unit,*) pcm%tess(ia)%point(3)
 
  876      do ia = 1, pcm%n_tesserae
 
  877        write(cav_gamess_unit,*) pcm%tess(ia)%area
 
  880      do ia = 1, pcm%n_tesserae
 
  881        write(cav_gamess_unit,*) pcm%tess(ia)%r_sphere
 
  884      do ia = 1, pcm%n_tesserae
 
  885        write(cav_gamess_unit,*) pcm%tess(ia)%normal
 
  893      safe_allocate(
mat_gamess(1:pcm%n_tesserae, 1:pcm%n_tesserae))
 
  897    if (.not. 
is_close(pcm%epsilon_infty, pcm%epsilon_0)) 
then 
  899      safe_allocate(pcm%matrix_d(1:pcm%n_tesserae, 1:pcm%n_tesserae))
 
  902      call pcm_matrix(pcm%epsilon_infty, pcm%tess, pcm%n_tesserae, pcm%matrix_d)
 
  905        pcmmat_gamess_unit = 
io_open(
pcm_dir//
'pcm_matrix_gamess_dyn.out', pcm%namespace, action=
'write')
 
  907        do jtess = 1, pcm%n_tesserae
 
  908          do itess = 1, pcm%n_tesserae
 
  909            write(pcmmat_gamess_unit,*) 
mat_gamess(itess,jtess) 
 
  919        safe_allocate(pcm%matrix_lf_d(1:pcm%n_tesserae, 1:pcm%n_tesserae))
 
  922        call pcm_matrix(pcm%epsilon_infty, pcm%tess, pcm%n_tesserae, pcm%matrix_lf_d, .
true.)
 
  926          pcmmat_unit = 
io_open(
pcm_dir//
'pcm_matrix_dynamic_lf.out', pcm%namespace, action=
'write')
 
  927          do jtess = 1, pcm%n_tesserae
 
  928            do itess = 1, pcm%n_tesserae
 
  929              write(pcmmat_unit,*) pcm%matrix_lf_d(itess,jtess)
 
  939    safe_allocate(pcm%matrix(1:pcm%n_tesserae, 1:pcm%n_tesserae))
 
  942    call pcm_matrix(pcm%epsilon_0, pcm%tess, pcm%n_tesserae, pcm%matrix)
 
  945      pcmmat_unit = 
io_open(
pcm_dir//
'pcm_matrix.out', pcm%namespace, action=
'write')
 
  947        pcm%namespace, action=
'write')
 
  949      do jtess = 1, pcm%n_tesserae
 
  950        do itess = 1, pcm%n_tesserae
 
  951          write(pcmmat_unit,*) pcm%matrix(itess,jtess)
 
  966      safe_allocate(pcm%matrix_lf(1:pcm%n_tesserae, 1:pcm%n_tesserae))
 
  969      call pcm_matrix(pcm%epsilon_0, pcm%tess, pcm%n_tesserae, pcm%matrix_lf, .
true.)
 
  973        pcmmat_unit = 
io_open(
pcm_dir//
'pcm_matrix_static_lf.out', pcm%namespace, action=
'write')
 
  974        do jtess = 1, pcm%n_tesserae
 
  975          do itess = 1, pcm%n_tesserae
 
  976            write(pcmmat_unit,*) pcm%matrix_lf(itess,jtess)
 
  984    message(1) = 
"Info: PCM response matrices has been evaluated" 
  999    call parse_variable(namespace, 
'PCMCalcMethod', pcm_calc_direct, pcm%calc_method)
 
 1006      do ia = 1, pcm%n_tesserae
 
 1007        if (pcm%tess(ia)%area > max_area) max_area = pcm%tess(ia)%area
 
 1011      default_nn = int(max_area*pcm%gaussian_width/minval(grid%spacing(1:pcm%space%dim)))
 
 1013      changed_default_nn = .false.
 
 1015      do ii = default_nn, 1, -1
 
 1020          changed_default_nn = .
true.
 
 1023      if (changed_default_nn) 
then 
 1024        call messages_write(
'PCM nearest neighbors have been reduced from ')
 
 1034      default_nn = pcm%tess_nn
 
 1048      call parse_variable(namespace, 
'PCMChargeSmearNN', default_nn, pcm%tess_nn)
 
 1058      safe_allocate(pcm%rho_n(1:grid%np_part))
 
 1059      safe_allocate(pcm%rho_e(1:grid%np_part))
 
 1060      if (pcm%localf) 
then 
 1061        safe_allocate(pcm%rho_ext(1:grid%np_part))
 
 1062        if (pcm%kick_is_present) 
then 
 1063          safe_allocate(pcm%rho_kick(1:grid%np_part))
 
 1069    safe_allocate(pcm%v_n(1:pcm%n_tesserae))
 
 1070    safe_allocate(pcm%q_n(1:pcm%n_tesserae))
 
 1071    safe_allocate(pcm%v_n_rs(1:grid%np))
 
 1076    safe_allocate(pcm%v_e(1:pcm%n_tesserae))
 
 1077    safe_allocate(pcm%q_e(1:pcm%n_tesserae))
 
 1078    safe_allocate(pcm%v_e_rs(1:grid%np))
 
 1083    safe_allocate(pcm%q_e_in(1:pcm%n_tesserae))
 
 1087    if (pcm%localf) 
then 
 1088      safe_allocate(pcm%v_ext(1:pcm%n_tesserae))
 
 1089      safe_allocate(pcm%q_ext(1:pcm%n_tesserae))
 
 1090      safe_allocate(pcm%v_ext_rs(1:grid%np))
 
 1095      safe_allocate(pcm%q_ext_in(1:pcm%n_tesserae))
 
 1098      if (pcm%kick_is_present) 
then 
 1099        safe_allocate(pcm%v_kick(1:pcm%n_tesserae))
 
 1100        safe_allocate(pcm%q_kick(1:pcm%n_tesserae))
 
 1101        safe_allocate(pcm%v_kick_rs(1:grid%np))
 
 1109    pcm%q_e_nominal = qtot
 
 1110    pcm%q_n_nominal = val_charge
 
 1121  subroutine pcm_calc_pot_rs(pcm, mesh, psolver, ions, v_h, v_ext, kick, time_present, kick_time)
 
 1123    type(
pcm_t),                
intent(inout) :: pcm
 
 1124    class(
mesh_t),              
intent(in)    :: mesh
 
 1126    type(
ions_t),     
optional, 
intent(in)    :: ions
 
 1127    real(real64),     
optional, 
intent(in)    :: v_h(:)
 
 1128    real(real64),     
optional, 
intent(in)    :: v_ext(:)
 
 1129    real(real64),     
optional, 
intent(in)    :: kick(:)
 
 1130    logical,          
optional, 
intent(in)    :: time_present
 
 1131    logical,          
optional, 
intent(in)    :: kick_time
 
 1135    logical :: input_asc_e
 
 1136    logical :: input_asc_ext
 
 1141    logical :: not_yet_called = .false.
 
 1142    logical :: is_time_for_kick = .false.
 
 1143    logical :: after_kick = .false.
 
 1145    logical :: td_calc_mode = .false.
 
 1147    integer :: asc_vs_t_unit_e
 
 1149    character(len=23) :: asc_vs_t_unit_format
 
 1150    character(len=16) :: asc_vs_t_unit_format_tail
 
 1156    assert(
present(v_h) .or. 
present(ions) .or. 
present(v_ext) .or. 
present(kick))
 
 1161    if ((.not. 
present(v_ext)) .and. 
present(kick)) calc = 
pcm_kick 
 1164    if (
present(time_present)) 
then 
 1165      if (time_present) td_calc_mode = .
true.
 
 1168    if (
present(kick_time)) 
then 
 1169      is_time_for_kick = kick_time
 
 1170      if (kick_time) after_kick = .
true.
 
 1173    select case (pcm%initial_asc)
 
 1175      input_asc_e = .
true.
 
 1176      input_asc_ext = .false.
 
 1178      input_asc_e = .false.
 
 1179      input_asc_ext = .
true.
 
 1181      input_asc_e = .
true.
 
 1182      input_asc_ext = .
true.
 
 1184      input_asc_e = .false.
 
 1185      input_asc_ext = .false.
 
 1190      call pcm_charges(pcm%q_n, pcm%qtot_n, pcm%v_n, pcm%matrix, pcm%n_tesserae, &
 
 1191        pcm%q_n_nominal, pcm%epsilon_0, pcm%renorm_charges, pcm%q_tot_tol, pcm%deltaQ_n)
 
 1195      call pcm_pot_rs(pcm, pcm%v_n_rs, pcm%q_n, pcm%rho_n, mesh, psolver)
 
 1200      if (td_calc_mode .and. .not. 
is_close(pcm%epsilon_infty, pcm%epsilon_0) .and. pcm%tdlevel /= 
pcm_td_eq) 
then 
 1202        select case (pcm%tdlevel)
 
 1204          select case (pcm%which_eps)
 
 1214          pcm%qtot_e = sum(pcm%q_e)
 
 1216          select case (pcm%iter)
 
 1221            call pcm_charges(pcm%q_e, pcm%qtot_e, pcm%v_e, pcm%matrix, pcm%n_tesserae, &
 
 1222              pcm%q_e_nominal, pcm%epsilon_0, pcm%renorm_charges, pcm%q_tot_tol, pcm%deltaQ_e)
 
 1225            call pcm_charges(pcm%q_e_in, pcm%qtot_e_in, pcm%v_e, pcm%matrix_d, pcm%n_tesserae, &
 
 1226              pcm%q_e_nominal, pcm%epsilon_infty, pcm%renorm_charges, pcm%q_tot_tol, pcm%deltaQ_e)
 
 1228            pcm%q_e_in = pcm%q_e - pcm%q_e_in
 
 1229            pcm%qtot_e_in = pcm%qtot_e - pcm%qtot_e_in
 
 1232            call pcm_charges(pcm%q_e, pcm%qtot_e, pcm%v_e, pcm%matrix_d, pcm%n_tesserae, &
 
 1233              pcm%q_e_nominal, pcm%epsilon_infty, pcm%renorm_charges, pcm%q_tot_tol, pcm%deltaQ_e)
 
 1235            pcm%q_e = pcm%q_e + pcm%q_e_in
 
 1236            pcm%qtot_e = pcm%qtot_e + pcm%qtot_e_in
 
 1242        call pcm_charges(pcm%q_e, pcm%qtot_e, pcm%v_e, pcm%matrix, pcm%n_tesserae, &
 
 1243          pcm%q_e_nominal, pcm%epsilon_0, pcm%renorm_charges, pcm%q_tot_tol, pcm%deltaQ_e)
 
 1249      call pcm_pot_rs(pcm, pcm%v_e_rs, pcm%q_e, pcm%rho_e, mesh, psolver)
 
 1259      if (td_calc_mode .and. .not. 
is_close(pcm%epsilon_infty, pcm%epsilon_0) .and. pcm%tdlevel /= 
pcm_td_eq) 
then 
 1261        select case (pcm%tdlevel)
 
 1263          select case (pcm%which_eps)
 
 1266              pcm%which_eps, pcm%namespace, this_drl=pcm%drl)
 
 1269              pcm%which_eps, pcm%namespace, this_deb=pcm%deb)
 
 1273          pcm%qtot_ext = sum(pcm%q_ext)
 
 1276          pcm%q_ext = pcm%q_ext + pcm%q_ext_in
 
 1277          pcm%qtot_ext = pcm%qtot_ext + pcm%qtot_ext_in
 
 1280          call pcm_charges(pcm%q_ext, pcm%qtot_ext, pcm%v_ext, pcm%matrix_lf_d, pcm%n_tesserae)
 
 1284          pcm%q_ext = pcm%q_ext + pcm%q_ext_in
 
 1285          pcm%qtot_ext = pcm%qtot_ext + pcm%qtot_ext_in
 
 1292        call pcm_charges(pcm%q_ext, pcm%qtot_ext, pcm%v_ext, pcm%matrix_lf, pcm%n_tesserae)
 
 1295        pcm%q_ext_in = pcm%q_ext
 
 1296        pcm%qtot_ext_in = pcm%qtot_ext
 
 1301      call pcm_pot_rs(pcm, pcm%v_ext_rs, pcm%q_ext, pcm%rho_ext, mesh, psolver)
 
 1308      if (is_time_for_kick) 
then 
 1313      if (pcm%kick_like) 
then 
 1314        if (is_time_for_kick) 
then 
 1316          if (.not. 
is_close(pcm%epsilon_infty, pcm%epsilon_0)) 
then 
 1317            call pcm_charges(pcm%q_kick, pcm%qtot_kick, pcm%v_kick, pcm%matrix_lf_d, pcm%n_tesserae)
 
 1319            call pcm_charges(pcm%q_kick, pcm%qtot_kick, pcm%v_kick, pcm%matrix_lf, pcm%n_tesserae)
 
 1331        if (after_kick) 
then 
 1333          select case (pcm%which_eps)
 
 1336              pcm%which_eps, pcm%namespace, this_drl=pcm%drl)
 
 1339              pcm%which_eps, pcm%namespace, this_deb=pcm%deb)
 
 1342          pcm%qtot_kick  = sum(pcm%q_kick)
 
 1358      call pcm_pot_rs(pcm, pcm%v_kick_rs, pcm%q_kick, pcm%rho_kick, mesh, psolver)
 
 1359      if (.not. pcm%kick_like) 
then 
 1361          pcm%q_ext    = pcm%q_ext    + pcm%q_kick
 
 1362          pcm%v_ext_rs = pcm%v_ext_rs + pcm%v_kick_rs
 
 1364          pcm%q_ext    = pcm%q_kick
 
 1365          pcm%v_ext_rs = pcm%v_kick_rs
 
 1373      write (asc_vs_t_unit_format_tail,
'(I5,A11)') pcm%n_tesserae,
'(1X,F14.8))' 
 1374      write (asc_vs_t_unit_format,
'(A)') 
'(F14.8,'//trim(adjustl(asc_vs_t_unit_format_tail))
 
 1376      if (pcm%solute .and. pcm%localf .and. td_calc_mode .and. calc == 
pcm_electrons) 
then 
 1377        asc_vs_t_unit_e = 
io_open(
pcm_dir//
'ASC_e_vs_t.dat', pcm%namespace, &
 
 1378          action=
'write', position=
'append', form=
'formatted')
 
 1379        write(asc_vs_t_unit_e,trim(adjustl(asc_vs_t_unit_format))) pcm%iter*pcm%dt, &
 
 1380          (pcm%q_e(ia) , ia = 1,pcm%n_tesserae)
 
 1392    type(
pcm_t), 
intent(in)  :: pcm
 
 1393    type(
mesh_t), 
intent(in) :: mesh
 
 1394    real(real64), 
intent(in)        :: v_mesh(:)
 
 1395    real(real64), 
intent(out)       :: v_cav(:)
 
 1407    do ia = 1, pcm%n_tesserae
 
 1421    real(real64),        
intent(out) :: v_n_cav(:)
 
 1422    type(
ions_t),        
intent(in)  :: ions
 
 1424    integer,             
intent(in)  :: n_tess
 
 1426    real(real64)   :: dist
 
 1434      do ia = 1, ions%natoms
 
 1436        dist = norm2(ions%pos(1:pcm_dim_space, ia) - tess(ik)%point)
 
 1438        v_n_cav(ik) = v_n_cav(ik) - ions%charge(ia)/dist
 
 1450  subroutine pcm_elect_energy(ions, pcm, E_int_ee, E_int_en, E_int_ne, E_int_nn, E_int_e_ext, E_int_n_ext)
 
 1451    type(
ions_t),     
intent(in)  :: ions
 
 1452    type(
pcm_t),      
intent(in)  :: pcm
 
 1453    real(real64),     
intent(out) :: e_int_ee
 
 1454    real(real64),     
intent(out) :: e_int_en
 
 1455    real(real64),     
intent(out) :: e_int_ne
 
 1456    real(real64),     
intent(out) :: e_int_nn
 
 1457    real(real64), 
optional,  
intent(out) :: e_int_e_ext
 
 1458    real(real64), 
optional,  
intent(out) :: e_int_n_ext
 
 1460    real(real64)   :: dist, z_ia
 
 1470    if (pcm%localf .and. ((.not. 
present(e_int_e_ext)) .or. &
 
 1471      (.not. 
present(e_int_n_ext)))) 
then 
 1472      message(1) = 
"pcm_elect_energy: There are lacking terms in subroutine call." 
 1474    else if (pcm%localf .and. (
present(e_int_e_ext) .and. &
 
 1475      present(e_int_n_ext))) 
then 
 1480    do ik = 1, pcm%n_tesserae
 
 1482      e_int_ee = e_int_ee + pcm%v_e(ik)*pcm%q_e(ik)
 
 1483      e_int_en = e_int_en + pcm%v_e(ik)*pcm%q_n(ik)
 
 1484      if (pcm%localf) 
then 
 1485        e_int_e_ext = e_int_e_ext + pcm%v_e(ik)*pcm%q_ext(ik)
 
 1488      do ia = 1, ions%natoms
 
 1490        dist = norm2(ions%pos(1:pcm_dim_space, ia) - pcm%tess(ik)%point)
 
 1492        z_ia = -ions%charge(ia)
 
 1494        e_int_ne = e_int_ne + z_ia*pcm%q_e(ik) / dist
 
 1495        e_int_nn = e_int_nn + z_ia*pcm%q_n(ik) / dist
 
 1496        if (pcm%localf) e_int_n_ext = e_int_n_ext + z_ia*pcm%q_ext(ik) / dist
 
 1501    e_int_ee = 
m_half*e_int_ee
 
 1502    e_int_en = 
m_half*e_int_en
 
 1503    e_int_ne = 
m_half*e_int_ne
 
 1504    e_int_nn = 
m_half*e_int_nn
 
 1509      if (pcm%localf) 
then 
 1510        write(pcm%info_unit, &
 
 1511          '(3X,I5,5X,F20.8,5X,F20.8,5X,F20.8,5X,F20.8,5X,F20.8,5X,F20.8,5X,F20.8,5X,F20.8,F20.8,5X,F20.8,5X,F20.8,5X,F20.8,5X)') &
 
 1526        write(pcm%info_unit, &
 
 1527          '(3X,I5,5X,F20.8,5X,F20.8,5X,F20.8,5X,F20.8,5X,F20.8,5X,F20.8,5X,F20.8,5X,F20.8,5X,F20.8)') &
 
 1549  subroutine pcm_charges(q_pcm, q_pcm_tot, v_cav, pcm_mat, n_tess, qtot_nominal, epsilon, renorm_charges, q_tot_tol, deltaQ)
 
 1550    real(real64),      
intent(out) :: q_pcm(:)
 
 1551    real(real64),      
intent(out) :: q_pcm_tot
 
 1552    real(real64),      
intent(in)  :: v_cav(:)
 
 1553    real(real64),      
intent(in)  :: pcm_mat(:,:)
 
 1554    integer,           
intent(in)  :: n_tess
 
 1555    real(real64),   
optional, 
intent(in)  :: qtot_nominal
 
 1556    real(real64),   
optional, 
intent(in)  :: epsilon
 
 1557    logical, 
optional, 
intent(in)  :: renorm_charges
 
 1558    real(real64),   
optional, 
intent(in)  :: q_tot_tol
 
 1559    real(real64),   
optional, 
intent(out) :: deltaq
 
 1562    real(real64) :: q_pcm_tot_norm
 
 1563    real(real64) :: coeff
 
 1573        q_pcm(ia) = q_pcm(ia) + pcm_mat(ia,ib)*v_cav(ib) 
 
 1575      q_pcm_tot = q_pcm_tot + q_pcm(ia)
 
 1578    if (
present(qtot_nominal)   .and. 
present(epsilon) .and.  &
 
 1579      present(renorm_charges) .and. 
present(q_tot_tol) .and. 
present(deltaq)) 
then 
 1581      deltaq = abs(q_pcm_tot) - ((epsilon - 
m_one)/epsilon) * abs(qtot_nominal)
 
 1582      if ((renorm_charges) .and. (abs(deltaq) > q_tot_tol)) 
then 
 1584        coeff = sign(
m_one, qtot_nominal)*sign(
m_one, deltaq)
 
 1586          q_pcm(ia) = q_pcm(ia) + coeff*q_pcm(ia)/q_pcm_tot*abs(deltaq)
 
 1587          q_pcm_tot_norm = q_pcm_tot_norm + q_pcm(ia)
 
 1589        q_pcm_tot = q_pcm_tot_norm
 
 1600    type(
pcm_t),   
intent(in) :: pcm
 
 1601    class(
mesh_t), 
intent(in) :: mesh
 
 1603    integer :: ia, nm(pcm%space%dim), ipt, i1, i2, i3
 
 1604    real(real64) :: posrel(pcm%space%dim)
 
 1605    integer(int64) :: pt
 
 1610    do ia = 1, pcm%n_tesserae
 
 1612      posrel(:) = pcm%tess(ia)%point(1:pcm%space%dim)/mesh%spacing(1:pcm%space%dim)
 
 1614      nm(:) = 
floor(posrel(:))
 
 1618      do i1 = -pcm%tess_nn + 1 , pcm%tess_nn
 
 1619        do i2 = -pcm%tess_nn + 1 , pcm%tess_nn
 
 1620          do i3 = -pcm%tess_nn + 1 , pcm%tess_nn
 
 1624            if (pt <= 0 .or. pt > mesh%np_part_global) 
then 
 1641    type(
pcm_t),      
intent(in) :: pcm
 
 1647      message(1) = 
'The simulation box is too small to contain all the requested' 
 1648      message(2) = 
'nearest neighbors for each tessera.' 
 1649      message(3) = 
'Consider using a larger box or reduce PCMChargeSmearNN.' 
 1660    type(
pcm_t),     
intent(inout) :: pcm
 
 1661    real(real64),    
intent(in)    :: q_pcm(:)
 
 1662    real(real64),    
intent(in)    :: q_pcm_tot
 
 1663    type(
mesh_t),    
intent(in)    :: mesh
 
 1664    real(real64),    
intent(out)   :: rho(:)
 
 1667    real(real64)   :: norm, qtot, rr, xx(pcm%space%dim), pp(pcm%space%dim)
 
 1670    integer :: nm(pcm%space%dim)
 
 1671    real(real64) :: posrel(pcm%space%dim)
 
 1673    integer :: i1, i2, i3
 
 1674    integer, 
allocatable :: pt(:)
 
 1675    real(real64),   
allocatable :: lrho(:) 
 
 1676    logical :: inner_point, boundary_point
 
 1684    npt = (2*pcm%tess_nn)**pcm%space%dim
 
 1685    safe_allocate(pt(1:npt))
 
 1686    safe_allocate(lrho(1:npt))
 
 1691    do ia = 1, pcm%n_tesserae
 
 1693      pp(:) = pcm%tess(ia)%point(1:pcm%space%dim)
 
 1694      posrel(:) = pp(:)/mesh%spacing(1:pcm%space%dim)
 
 1696      nm(:) = 
floor(posrel(:))
 
 1700      do i1 = -pcm%tess_nn + 1 , pcm%tess_nn
 
 1701        do i2 = -pcm%tess_nn + 1 , pcm%tess_nn
 
 1702          do i3 = -pcm%tess_nn + 1 , pcm%tess_nn
 
 1718        if (pt(ipt) > 0 .and. pt(ipt) <= mesh%np_part) 
then 
 1720          if (mesh%parallel_in_domains) 
then 
 1721            boundary_point = pt(ipt) > mesh%np + mesh%pv%np_ghost
 
 1722            inner_point = pt(ipt) > 0 .and. pt(ipt) <= mesh%np
 
 1724            if (boundary_point .or. inner_point) 
then 
 1725              xx(:) = mesh%x(pt(ipt),:)
 
 1731            xx(:) = mesh%x(pt(ipt), :)
 
 1734          rr = sum((xx(1:pcm%space%dim) - pp(1:pcm%space%dim))**2)
 
 1735          norm = norm + 
exp(-rr/(pcm%tess(ia)%area*pcm%gaussian_width))
 
 1736          lrho(ipt) = lrho(ipt) + 
exp(-rr/(pcm%tess(ia)%area*pcm%gaussian_width))
 
 1742      call mesh%allreduce(lrho, npt)
 
 1745      norm = sum(lrho(1:npt)) * mesh%volume_element
 
 1747        norm = q_pcm(ia)/norm
 
 1751      lrho(:) = lrho(:) * norm
 
 1756        if (pt(ipt) > 0 .and. pt(ipt) <= mesh%np_part_global) 
then 
 1758          if (mesh%parallel_in_domains) 
then 
 1759            boundary_point = pt(ipt) > mesh%np + mesh%pv%np_ghost
 
 1760            inner_point = pt(ipt) > 0 .and. pt(ipt) <= mesh%np
 
 1761            if (boundary_point .or. inner_point) rho(pt(ipt)) = rho(pt(ipt)) + lrho(ipt)
 
 1763            rho(pt(ipt)) = rho(pt(ipt)) + lrho(ipt)
 
 1771    if (
debug%info) 
then 
 1785    safe_deallocate_a(pt)
 
 1786    safe_deallocate_a(lrho)
 
 1796  subroutine pcm_pot_rs(pcm, v_pcm, q_pcm, rho, mesh, psolver)
 
 1797    type(
pcm_t),       
intent(inout) :: pcm
 
 1798    real(real64), 
contiguous, 
intent(inout) :: v_pcm(:)
 
 1799    real(real64), 
contiguous, 
intent(in)    :: q_pcm(:)
 
 1800    real(real64), 
contiguous, 
intent(inout) :: rho(:)
 
 1801    type(
mesh_t),      
intent(in)    :: mesh
 
 1811    select case (pcm%calc_method)
 
 1812    case (pcm_calc_direct)
 
 1813      call pcm_pot_rs_direct(v_pcm, q_pcm, pcm%tess, pcm%n_tesserae, mesh, pcm%gaussian_width)
 
 1824    if (
debug%info) 
then 
 1839    real(real64), 
contiguous, 
intent(inout) :: v_pcm(:)
 
 1841    real(real64), 
contiguous, 
intent(inout) :: rho(:)
 
 1854    real(real64),        
intent(out) :: v_pcm(:)
 
 1855    real(real64),        
intent(in)  :: q_pcm(:)
 
 1856    real(real64),        
intent(in)  :: width_factor
 
 1857    integer,             
intent(in)  :: n_tess
 
 1858    type(
mesh_t),        
intent(in)  :: mesh
 
 1861    real(real64), 
parameter :: p_1 = 0.119763_real64
 
 1862    real(real64), 
parameter :: p_2 = 0.205117_real64
 
 1863    real(real64), 
parameter :: q_1 = 0.137546_real64
 
 1864    real(real64), 
parameter :: q_2 = 0.434344_real64
 
 1866    real(real64)     :: term
 
 1879          call mesh_r(mesh, ip, term, origin=tess(ia)%point)
 
 1880          arg = term/
sqrt(tess(ia)%area*width_factor)
 
 1881          term = (
m_one + p_1*arg + p_2*arg**2) / (
m_one + q_1*arg + q_2*arg**2 + p_2*arg**3)
 
 1882          v_pcm(ip) = v_pcm(ip) + q_pcm(ia) * term/
sqrt(tess(ia)%area*width_factor)
 
 1893          call mesh_r(mesh, ip, term, origin=tess(ia)%point)
 
 1894          v_pcm(ip) = v_pcm(ip) + q_pcm(ia)/term
 
 1905  subroutine pcm_matrix(eps, tess, n_tess, pcm_mat, localf)
 
 1906    real(real64),        
intent(in)  :: eps
 
 1908    integer,             
intent(in)  :: n_tess
 
 1909    real(real64),        
intent(out) :: pcm_mat(:,:)
 
 1910    logical,   
optional, 
intent(in)  :: localf
 
 1913    integer, 
allocatable :: iwork(:)
 
 1914    real(real64), 
allocatable :: mat_tmp(:,:)
 
 1916    real(real64) :: sgn_lf
 
 1921    safe_allocate(s_mat_act(1:n_tess, 1:n_tess))
 
 1925    safe_allocate(sigma(1:n_tess, 1:n_tess))
 
 1926    sigma = s_mat_act/eps
 
 1929    safe_allocate(d_mat_act(1:n_tess, 1:n_tess))
 
 1933    safe_allocate(delta(1:n_tess, 1:n_tess))
 
 1938    if (
present(localf)) 
then 
 1939      if (localf) sgn_lf = -
m_one 
 1943    pcm_mat = - sgn_lf * d_mat_act
 
 1949    safe_deallocate_a(d_mat_act)
 
 1951    safe_allocate(iwork(1:n_tess))
 
 1956    call dgesv(n_tess, n_tess, s_mat_act, n_tess, iwork, pcm_mat, n_tess, 
info)
 
 1958    safe_deallocate_a(iwork)
 
 1960    safe_deallocate_a(s_mat_act)
 
 1964    pcm_mat = -matmul(sigma, pcm_mat)
 
 1967      pcm_mat(i,i) = pcm_mat(i,i) + 
m_two*
m_pi 
 1970    pcm_mat = pcm_mat - sgn_lf * delta
 
 1972    safe_allocate(mat_tmp(1:n_tess,1:n_tess))
 
 1975    safe_allocate(d_mat_act(1:n_tess,1:n_tess))
 
 1978    mat_tmp = transpose(d_mat_act)
 
 1980    mat_tmp = matmul(sigma, mat_tmp)
 
 1984    safe_deallocate_a(d_mat_act)
 
 1986    safe_allocate(s_mat_act(1:n_tess, 1:n_tess))
 
 1989    mat_tmp = mat_tmp + 
m_two*
m_pi*s_mat_act - matmul(delta, s_mat_act)
 
 1991    safe_deallocate_a(s_mat_act)
 
 1992    safe_deallocate_a(sigma)
 
 1993    safe_deallocate_a(delta)
 
 1995    safe_allocate(iwork(1:n_tess))
 
 1999    call dgesv(n_tess, n_tess, mat_tmp, n_tess, iwork, pcm_mat, n_tess, 
info)
 
 2001    safe_deallocate_a(iwork)
 
 2002    safe_deallocate_a(mat_tmp)
 
 2004    pcm_mat = - sgn_lf * pcm_mat
 
 2019    integer,             
intent(in)    :: n_tess
 
 2029        if (ii /= jj) s_mat_act(jj,ii) = s_mat_act(ii,jj) 
 
 2038    integer,             
intent(in) :: n_tess
 
 2057  real(real64) function s_mat_elem_I(tessi, tessj)
 
 2061    real(real64), 
parameter :: M_SD_DIAG    = 1.0694_real64
 
 2062    real(real64), 
parameter :: M_DIST_MIN   = 0.1_real64
 
 2064    real(real64) :: diff(1:PCM_DIM_SPACE)
 
 2065    real(real64) :: dist
 
 2066    real(real64) :: s_diag
 
 2067    real(real64) :: s_off_diag
 
 2072    diff = tessi%point - tessj%point
 
 2078      s_mat_elem_i = s_diag
 
 2080      if (dist > m_dist_min) s_off_diag = 
m_one/dist
 
 2081      s_mat_elem_i = s_off_diag
 
 2089  real(real64) function d_mat_elem_I(tessi, tessj)
 
 2090    type(pcm_tessera_t), 
intent(in) :: tessi
 
 2091    type(pcm_tessera_t), 
intent(in) :: tessj
 
 2093    real(real64), 
parameter :: M_SD_DIAG    = 1.0694_real64
 
 2094    real(real64), 
parameter :: M_DIST_MIN   = 0.04_real64
 
 2096    real(real64) :: diff(1:PCM_DIM_SPACE)
 
 2097    real(real64) :: dist
 
 2098    real(real64) :: d_diag
 
 2099    real(real64) :: d_off_diag
 
 2105    diff = tessi%point - tessj%point
 
 2109    if (abs(dist) <= m_epsilon) 
then 
 2111      d_diag = m_sd_diag*
sqrt(m_four*m_pi*tessi%area)
 
 2112      d_diag = -d_diag/(m_two*tessi%r_sphere)
 
 2117      if (dist > m_dist_min) 
then 
 2118        d_off_diag = dot_product( diff, tessj%normal(:))
 
 2119        d_off_diag = d_off_diag*tessj%area/dist**3
 
 2133  subroutine cav_gen(tess_sphere, tess_min_distance, nesf, sfe, nts, cts, unit_pcminfo)
 
 2134    integer,              
intent(in)    :: tess_sphere
 
 2135    real(real64)  ,              
intent(in)    :: tess_min_distance
 
 2137    integer,              
intent(in)    :: nesf
 
 2138    integer,              
intent(out)   :: nts
 
 2139    type(pcm_tessera_t),  
intent(out)   :: cts(:)
 
 2140    integer,              
intent(in)    :: unit_pcminfo
 
 2142    integer, 
parameter :: DIM_ANGLES = 24
 
 2143    integer, 
parameter :: DIM_TEN = 10
 
 2144    integer, 
parameter :: DIM_VERTICES = 122
 
 2145    integer, 
parameter :: MAX_VERTICES = 6
 
 2146    integer, 
parameter :: MXTS = 10000
 
 2148    real(real64) :: thev(1:DIM_ANGLES)
 
 2149    real(real64) :: fiv(1:DIM_ANGLES)
 
 2170    integer :: isfet(1:dim_ten*dim_angles)
 
 2188    real(real64) :: area
 
 2189    real(real64) :: test2
 
 2191    real(real64) :: dnorm
 
 2201    real(real64) :: stot
 
 2202    real(real64) :: prod
 
 2205    logical :: band_iter
 
 2211    data thev/ 0.6523581398_real64 , 1.107148718_real64  , 1.382085796_real64 , &
 
 2212      1.759506858_real64  , 2.034443936_real64  , 2.489234514_real64 , &
 
 2213      0.3261790699_real64 , 0.5535743589_real64, &
 
 2214      0.8559571251_real64 , 0.8559571251_real64 , 1.017221968_real64 , &
 
 2215      1.229116717_real64  , 1.229116717_real64  , 1.433327788_real64 , &
 
 2216      1.570796327_real64  , 1.570796327_real64  , 1.708264866_real64 , &
 
 2217      1.912475937_real64  , 1.912475937_real64  , 2.124370686_real64 , &
 
 2218      2.285635528_real64  , 2.285635528_real64  , 2.588018295_real64 , &
 
 2219      2.815413584_real64 /
 
 2220    data fiv/                       0.6283185307_real64 , m_zero            , &
 
 2221      0.6283185307_real64 , m_zero             , 0.6283185307_real64, &
 
 2222      m_zero             , 0.6283185307_real64 , m_zero,       &
 
 2223      0.2520539002_real64 , 1.004583161_real64  , 0.6283185307_real64, &
 
 2224      0.3293628477_real64 , 0.9272742138_real64 , m_zero,        &
 
 2225      0.3141592654_real64 , 0.9424777961_real64 , 0.6283185307_real64, &
 
 2226      0.2989556830_real64 , 0.9576813784_real64 , m_zero,        &
 
 2227      0.3762646305_real64 , 0.8803724309_real64 , 0.6283188307_real64, &
 
 2229    data fir / 1.256637061_real64  /
 
 2233    data (idum(ii),ii = 1, 280) /                                   &
 
 2234      1, 6, 2, 32, 36, 37, 1, 2, 3, 33, 32, 38, 1, 3, 4, 34,         &
 
 2235      33, 39, 1, 4, 5, 35, 34, 40, 1, 5, 6, 36, 35, 41, 7, 2, 6, 51, &
 
 2236      42, 37, 8, 3, 2, 47, 43, 38, 9, 4, 3, 48, 44, 39, 10, 5, 4,    &
 
 2237      49, 45, 40, 11, 6, 5, 50, 46, 41, 8, 2, 12, 62, 47, 52, 9,     &
 
 2238      3, 13, 63, 48, 53, 10, 4, 14, 64, 49, 54, 11, 5, 15, 65, 50,   &
 
 2239      55, 7, 6, 16, 66, 51, 56, 7, 12, 2, 42, 57, 52, 8, 13, 3,      &
 
 2240      43, 58, 53, 9, 14, 4, 44, 59, 54, 10, 15, 5, 45, 60, 55, 11,   &
 
 2241      16, 6, 46, 61, 56, 8, 12, 18, 68, 62, 77, 9, 13, 19, 69, 63,   &
 
 2242      78, 10, 14, 20, 70, 64, 79, 11, 15, 21, 71, 65, 80, 7, 16,     &
 
 2243      17, 67, 66, 81, 7, 17, 12, 57, 67, 72, 8, 18, 13, 58, 68, 73,  &
 
 2244      9, 19, 14, 59, 69, 74, 10, 20, 15, 60, 70, 75, 11, 21, 16,     &
 
 2245      61, 71, 76, 22, 12, 17, 87, 82, 72, 23, 13, 18, 88, 83, 73,    &
 
 2246      24, 14, 19, 89, 84, 74, 25, 15, 20, 90, 85, 75, 26, 16, 21,    &
 
 2247      91, 86, 76, 22, 18, 12, 82, 92, 77, 23, 19, 13, 83, 93, 78,    &
 
 2248      24, 20, 14, 84, 94, 79, 25, 21, 15, 85, 95, 80, 26, 17, 16,    &
 
 2249      86, 96, 81, 22, 17, 27, 102, 87, 97, 23, 18, 28, 103, 88, 98,  &
 
 2250      24, 19, 29, 104, 89, 99, 25, 20, 30, 105, 90, 100, 26, 21,     &
 
 2251      31, 106, 91, 101, 22, 28, 18, 92, 107, 98, 23, 29, 19, 93 /
 
 2252    data (idum(ii),ii = 281,360) /              &
 
 2253      108, 99, 24, 30, 20, 94, 109, 100, 25, 31, 21, 95, 110, 101,   &
 
 2254      26, 27, 17, 96, 111, 97, 22, 27, 28, 107, 102, 112, 23, 28,    &
 
 2255      29, 108, 103, 113, 24, 29, 30, 109, 104, 114, 25, 30, 31,      &
 
 2256      110, 105, 115, 26, 31, 27, 111, 106, 116, 122, 28, 27, 117,    &
 
 2257      118, 112, 122, 29, 28, 118, 119, 113, 122, 30, 29, 119, 120,   &
 
 2258      114, 122, 31, 30, 120, 121, 115, 122, 27, 31, 121, 117, 116 /
 
 2260    if (mpi_grp_is_root(mpi_world)) 
then 
 2261      if (tess_sphere == 1) 
then 
 2262        write(unit_pcminfo,
'(A1)')  
'#' 
 2263        write(unit_pcminfo,
'(A34)') 
'# Number of tesserae / sphere = 60' 
 2264        write(unit_pcminfo,
'(A1)')  
'#' 
 2266        write(unit_pcminfo,
'(A1)')  
'#' 
 2267        write(unit_pcminfo,
'(A35)') 
'# Number of tesserae / sphere = 240' 
 2268        write(unit_pcminfo,
'(A1)')  
'#' 
 2277    sfe(:)%x = sfe(:)%x*p_a_b
 
 2278    sfe(:)%y = sfe(:)%y*p_a_b
 
 2279    sfe(:)%z = sfe(:)%z*p_a_b
 
 2280    sfe(:)%r = sfe(:)%r*p_a_b
 
 2284    jvt1 = reshape(idum,(/6,60/))
 
 2307    do ia = 1, dim_angles
 
 2314        if (ja == 1) fi = fiv(ia)
 
 2316        cv(ii,1) = sth*
cos(fi)
 
 2317        cv(ii,2) = sth*
sin(fi)
 
 2336        do i_tes = 1, tess_sphere
 
 2337          if (tess_sphere == 1) 
then 
 2342            if (i_tes == 1)     
then 
 2346            elseif (i_tes == 2) 
then 
 2350            elseif (i_tes == 3) 
then 
 2354            elseif (i_tes == 4) 
then 
 2361          pts(1,1) = cv(n1,1)*ren + xen
 
 2362          pts(2,1) = cv(n1,3)*ren + yen
 
 2363          pts(3,1) = cv(n1,2)*ren + zen
 
 2365          pts(1,2) = cv(n2,1)*ren + xen
 
 2366          pts(2,2) = cv(n2,3)*ren + yen
 
 2367          pts(3,2) = cv(n2,2)*ren + zen
 
 2369          pts(1,3) = cv(n3,1)*ren + xen
 
 2370          pts(2,3) = cv(n3,3)*ren + yen
 
 2371          pts(3,3) = cv(n3,2)*ren + zen
 
 2377          call subtessera(sfe, nsfe, nesf, nv, pts ,ccc, pp, pp1, area)
 
 2379          if (abs(area) <= m_epsilon) cycle
 
 2381          xctst(tess_sphere*(its-1) + i_tes)   = pp(1)
 
 2382          yctst(tess_sphere*(its-1) + i_tes)   = pp(2)
 
 2383          zctst(tess_sphere*(its-1) + i_tes)   = pp(3)
 
 2384          nctst(:,tess_sphere*(its-1) + i_tes) = pp1(:)
 
 2385          ast(tess_sphere*(its-1) + i_tes)     = area
 
 2386          isfet(tess_sphere*(its-1) + i_tes)   = nsfe
 
 2393        if (abs(ast(its)) <= m_epsilon) cycle
 
 2397          write(message(1),
'(a,I5,a,I5)') 
"total number of tesserae", nn, 
">",mxts
 
 2398          call messages_warning(1)
 
 2401        cts(nn)%point(1)  = xctst(its)
 
 2402        cts(nn)%point(2)  = yctst(its)
 
 2403        cts(nn)%point(3)  = zctst(its)
 
 2404        cts(nn)%normal(:) = nctst(:,its)
 
 2405        cts(nn)%area      = ast(its)
 
 2406        cts(nn)%r_sphere  = sfe(isfet(its))%r
 
 2414    test2 = tess_min_distance*tess_min_distance
 
 2417    do while (.not.(band_iter))
 
 2420      loop_ia: 
do ia = 1, nts-1
 
 2421        if (abs(cts(ia)%area) <= m_epsilon) cycle
 
 2422        xi = cts(ia)%point(1)
 
 2423        yi = cts(ia)%point(2)
 
 2424        zi = cts(ia)%point(3)
 
 2426        loop_ja: 
do ja = ia+1, nts
 
 2427          if (abs(cts(ja)%area) <= m_epsilon) cycle
 
 2428          xj = cts(ja)%point(1)
 
 2429          yj = cts(ja)%point(2)
 
 2430          zj = cts(ja)%point(3)
 
 2432          rij = (xi - xj)**2 + (yi - yj)**2 + (zi - zj)**2
 
 2434          if (rij > test2) cycle
 
 2436          if (mpi_grp_is_root(mpi_world)) 
then 
 2437            write(unit_pcminfo,
'(A40,I4,A5,I4,A4,F8.4,A13,F8.4,A3)') &
 
 2438              '# Warning: The distance between tesserae', &
 
 2439              ia,
' and ', ja,
' is ',
sqrt(rij),
' A, less than', tess_min_distance,
' A.' 
 2443          xi = (xi*cts(ia)%area + xj*cts(ja)%area) / (cts(ia)%area + cts(ja)%area)
 
 2444          yi = (yi*cts(ia)%area + yj*cts(ja)%area) / (cts(ia)%area + cts(ja)%area)
 
 2445          zi = (zi*cts(ia)%area + zj*cts(ja)%area) / (cts(ia)%area + cts(ja)%area)
 
 2447          cts(ia)%point(1) = xi
 
 2448          cts(ia)%point(2) = yi
 
 2449          cts(ia)%point(3) = zi
 
 2452          cts(ia)%normal = (cts(ia)%normal*cts(ia)%area + cts(ja)%normal*cts(ja)%area)
 
 2453          dnorm = norm2(cts(ia)%normal)
 
 2454          cts(ia)%normal = cts(ia)%normal/dnorm
 
 2457          cts(ia)%r_sphere = (cts(ia)%r_sphere*cts(ia)%area + cts(ja)%r_sphere*cts(ja)%area) / &
 
 2458            (cts(ia)%area + cts(ja)%area)
 
 2461          cts(ia)%area = cts(ia)%area + cts(ja)%area
 
 2478      prod = dot_product(cts(its)%point, cts(its)%normal)
 
 2479      vol  = vol + cts(its)%area * prod / m_three
 
 2480      stot = stot + cts(its)%area
 
 2483    if (mpi_grp_is_root(mpi_world)) 
then 
 2484      write(unit_pcminfo, 
'(A2)')  
'# ' 
 2485      write(unit_pcminfo, 
'(A29,I4)')    
'# Total number of tesserae = ' , nts
 
 2486      write(unit_pcminfo, 
'(A30,F12.6)') 
'# Cavity surface area (A^2) = ', stot
 
 2487      write(unit_pcminfo, 
'(A24,F12.6)') 
'# Cavity volume (A^3) = '      , vol
 
 2488      write(unit_pcminfo, 
'(A2)')  
'# ' 
 2492    cts(:)%area     = cts(:)%area*p_ang*p_ang
 
 2493    cts(:)%point(1) = cts(:)%point(1)*p_ang
 
 2494    cts(:)%point(2) = cts(:)%point(2)*p_ang
 
 2495    cts(:)%point(3) = cts(:)%point(3)*p_ang
 
 2496    cts(:)%r_sphere = cts(:)%r_sphere*p_ang
 
 2498    sfe(:)%x=sfe(:)%x*p_ang
 
 2499    sfe(:)%y=sfe(:)%y*p_ang
 
 2500    sfe(:)%z=sfe(:)%z*p_ang
 
 2501    sfe(:)%r=sfe(:)%r*p_ang
 
 2510  subroutine subtessera(sfe, ns, nesf, nv, pts, ccc, pp, pp1, area)
 
 2512    integer,            
intent(in)    :: ns
 
 2513    integer,            
intent(in)    :: nesf
 
 2514    integer,            
intent(inout) :: nv
 
 2515    real(real64),       
intent(inout) :: pts(:,:)
 
 2516    real(real64),       
intent(out)   :: ccc(:,:)
 
 2517    real(real64),       
intent(out)   :: pp(:)
 
 2518    real(real64),       
intent(out)   :: pp1(:)
 
 2519    real(real64),       
intent(out)   :: area
 
 2521    real(real64), 
parameter   :: TOL = -1.0e-10_real64
 
 2522    integer, 
parameter :: DIM_TEN = 10
 
 2524    integer :: intsph(1:DIM_TEN)
 
 2535    real(real64)  :: p1(1:PCM_DIM_SPACE)
 
 2536    real(real64)  :: p2(1:PCM_DIM_SPACE)
 
 2537    real(real64)  :: p3(1:PCM_DIM_SPACE)
 
 2538    real(real64)  :: p4(1:PCM_DIM_SPACE)
 
 2539    real(real64)  :: point(1:PCM_DIM_SPACE)
 
 2540    real(real64)  :: pscr(1:PCM_DIM_SPACE,1:DIM_TEN)
 
 2541    real(real64)  :: cccp(1:PCM_DIM_SPACE,1:DIM_TEN)
 
 2542    real(real64)  :: pointl(1:PCM_DIM_SPACE,1:DIM_TEN)
 
 2543    real(real64)  :: diff(1:PCM_DIM_SPACE)
 
 2545    integer :: ind(1:DIM_TEN)
 
 2546    integer :: ltyp(1:DIM_TEN)
 
 2547    integer :: intscr(1:DIM_TEN)
 
 2549    real(real64)  :: delr
 
 2550    real(real64)  :: delr2
 
 2553    real(real64)  :: dnorm
 
 2554    real(real64)  :: dist
 
 2571      ccc(1,jj) = sfe(ns)%x
 
 2572      ccc(2,jj) = sfe(ns)%y
 
 2573      ccc(3,jj) = sfe(ns)%z
 
 2578      if (nsfe1 == ns) cycle
 
 2580        intscr(jj) = intsph(jj)
 
 2581        pscr(:,jj) = pts(:,jj)
 
 2582        cccp(:,jj) = ccc(:,jj)
 
 2590        delr2 = (pts(1,ii) - sfe(nsfe1)%x)**2 + (pts(2,ii) - sfe(nsfe1)%y)**2 + (pts(3,ii) - sfe(nsfe1)%z)**2
 
 2592        if (delr < sfe(nsfe1)%r) 
then 
 2598      if (icop == nv) 
then 
 2606        if (ll == nv) iv2 = 1
 
 2607        if ((ind(iv1) == 1) .and. (ind(iv2) == 1)) 
then 
 2609        else if ((ind(iv1) == 0) .and. (ind(iv2) == 1)) 
then 
 2611        else if ((ind(iv1) == 1) .and. (ind(iv2) == 0)) 
then 
 2613        else if ((ind(iv1) == 0) .and. (ind(iv2) == 0)) 
then 
 2615          diff = ccc(:,ll) - pts(:,ll)
 
 2616          rc2 = dot_product(diff, diff)
 
 2620            point = pts(:,iv1) + ii * (pts(:,iv2) - pts(:,iv1)) / 11
 
 2621            point = point - ccc(:,ll)
 
 2622            dnorm = norm2(point)
 
 2623            point = point * rc / dnorm + ccc(:,ll)
 
 2625            dist = 
sqrt((point(1) - sfe(nsfe1)%x)**2 + (point(2) - sfe(nsfe1)%y)**2 + (point(3) - sfe(nsfe1)%z)**2)
 
 2627            if ((dist - sfe(nsfe1)%r) < tol) 
then 
 2629              pointl(:, ll) = point
 
 2639        if ((ltyp(ll) == 1) .or. (ltyp(ll) == 2)) icut = icut + 1
 
 2640        if (ltyp(ll) == 3) icut = icut + 2
 
 2651        if (ltyp(ll) == 0) cycle
 
 2654        if (ll == nv) iv2 = 1
 
 2656        if (ltyp(ll) == 1) 
then 
 2657          pts(:,na) = pscr(:,iv1)
 
 2658          ccc(:,na) = cccp(:,iv1)
 
 2659          intsph(na) = intscr(iv1)
 
 2665          call inter(sfe, p1, p2, p3, p4, nsfe1, 0)
 
 2668          de2 = (sfe(nsfe1)%x - sfe(ns)%x)**2 + ( sfe(nsfe1)%y - sfe(ns)%y)**2 + &
 
 2669            (sfe(nsfe1)%z - sfe(ns)%z)**2
 
 2671          ccc(1,na) = sfe(ns)%x + ( sfe(nsfe1)%x - sfe(ns)%x)* &
 
 2672            (sfe(ns)%r**2 - sfe(nsfe1)%r**2 + de2) / (m_two*de2)
 
 2674          ccc(2,na) = sfe(ns)%y + ( sfe(nsfe1)%y - sfe(ns)%y)* &
 
 2675            (sfe(ns)%r**2 - sfe(nsfe1)%r**2 + de2) / (m_two*de2)
 
 2677          ccc(3,na) = sfe(ns)%z + ( sfe(nsfe1)%z - sfe(ns)%z)* &
 
 2678            (sfe(ns)%r**2 - sfe(nsfe1)%r**2 + de2) / (m_two*de2)
 
 2684        if (ltyp(ll) == 2) 
then 
 2689          call inter(sfe, p1, p2, p3, p4, nsfe1, 1)
 
 2691          ccc(:,na) = cccp(:,iv1)
 
 2692          intsph(na) = intscr(iv1)
 
 2696        if (ltyp(ll) == 3) 
then 
 2697          pts(:,na) = pscr(:,iv1)
 
 2698          ccc(:,na) = cccp(:,iv1)
 
 2699          intsph(na) = intscr(iv1)
 
 2705          call inter(sfe, p1, p2, p3, p4, nsfe1, 0)
 
 2708          de2 = (sfe(nsfe1)%x - sfe(ns)%x)**2 + (sfe(nsfe1)%y - sfe(ns)%y)**2 + (sfe(nsfe1)%z - sfe(ns)%z)**2
 
 2710          ccc(1,na) = sfe(ns)%x + (sfe(nsfe1)%x - sfe(ns)%x) * (sfe(ns)%r**2 - sfe(nsfe1)%r**2 + de2) / (m_two*de2)
 
 2712          ccc(2,na) = sfe(ns)%y + (sfe(nsfe1)%y - sfe(ns)%y) * (sfe(ns)%r**2 - sfe(nsfe1)%r**2 + de2) / (m_two*de2)
 
 2714          ccc(3,na) = sfe(ns)%z + (sfe(nsfe1)%z - sfe(ns)%z) * (sfe(ns)%r**2 - sfe(nsfe1)%r**2 + de2) / (m_two*de2)
 
 2722          call inter(sfe, p1, p2, p3, p4, nsfe1, 1)
 
 2724          ccc(:,na) = cccp(:,iv1)
 
 2725          intsph(na) = intscr(iv1)
 
 2729        if (ltyp(ll) == 4) 
then 
 2730          pts(:,na) = pscr(:,iv1)
 
 2731          ccc(:,na) = cccp(:,iv1)
 
 2732          intsph(na) = intscr(iv1)
 
 2739        message(1) = 
"Too many vertices on the tessera" 
 2740        call messages_fatal(1)
 
 2744    call gaubon(sfe, nv, ns, pts, ccc, pp, pp1, area, intsph)
 
 2754  subroutine inter(sfe, p1, p2, p3, p4, ns, ia)
 
 2756    real(real64),       
intent(in)  :: p1(1:PCM_DIM_SPACE)
 
 2757    real(real64),       
intent(in)  :: p2(1:PCM_DIM_SPACE)
 
 2758    real(real64),       
intent(in)  :: p3(1:PCM_DIM_SPACE)
 
 2759    real(real64),       
intent(out) :: p4(1:PCM_DIM_SPACE)
 
 2760    integer,            
intent(in)  :: ns
 
 2761    integer,            
intent(in)  :: ia
 
 2763    real(real64), 
parameter :: TOL = 1.0e-08_real64
 
 2767    real(real64)  :: alpha
 
 2768    real(real64)  :: delta
 
 2769    real(real64)  :: dnorm
 
 2770    real(real64)  :: diff
 
 2771    real(real64)  :: diff_vec(1:PCM_DIM_SPACE)
 
 2772    logical :: band_iter
 
 2784    do while(.not.(band_iter))
 
 2785      if (m_iter > 1000) 
then 
 2786        message(1) = 
"Too many iterations inside subrotuine inter" 
 2787        call messages_fatal(1)
 
 2792      alpha = alpha + delta
 
 2794      p4 = p1 + alpha*(p2 - p1) - p3
 
 2796      p4 = p4*r/dnorm + p3
 
 2797      diff = (p4(1) - sfe(ns)%x)**2 + (p4(2) - sfe(ns)%y)**2 + (p4(3) - sfe(ns)%z)**2
 
 2798      diff = 
sqrt(diff) - sfe(ns)%r
 
 2800      if (abs(diff) < tol) 
then 
 2806        if (diff > m_zero) delta =  m_one/(m_two**(m_iter + 1))
 
 2807        if (diff < m_zero) delta = -m_one/(m_two**(m_iter + 1))
 
 2813        if (diff > m_zero) delta = -m_one/(m_two**(m_iter + 1))
 
 2814        if (diff < m_zero) delta =  m_one/(m_two**(m_iter + 1))
 
 2821  end subroutine inter 
 2831  subroutine gaubon(sfe, nv, ns, pts, ccc, pp, pp1, area, intsph)
 
 2833    real(real64),       
intent(in)    :: pts(:,:)
 
 2834    real(real64),       
intent(in)    :: ccc(:,:)
 
 2835    real(real64),       
intent(inout) :: pp(:)
 
 2836    real(real64),       
intent(inout) :: pp1(:)
 
 2837    integer,            
intent(in)    :: intsph(:)
 
 2838    real(real64),       
intent(out)   :: area
 
 2839    integer,            
intent(in)    :: nv
 
 2840    integer,            
intent(in)    :: ns
 
 2842    real(real64) :: p1(1:PCM_DIM_SPACE), p2(1:PCM_DIM_SPACE), p3(1:PCM_DIM_SPACE)
 
 2843    real(real64) :: u1(1:PCM_DIM_SPACE), u2(1:PCM_DIM_SPACE)
 
 2844    real(real64) :: point_1(1:PCM_DIM_SPACE), point_2(1:PCM_DIM_SPACE)
 
 2845    real(real64) :: tpi, sum1, dnorm, dnorm1, dnorm2
 
 2846    real(real64) :: cosphin, phin, costn, sum2, betan
 
 2847    integer :: nsfe1, ia, nn, n0, n1, n2
 
 2862      point_1 = pts(:,nn) - ccc(:,nn)
 
 2864        point_2 = pts(:,nn+1) - ccc(:,nn)
 
 2866        point_2 = pts(:,1) - ccc(:,nn)
 
 2869      dnorm1 = norm2(point_1)
 
 2870      dnorm2 = norm2(point_2)
 
 2871      cosphin = dot_product(point_1, point_2) / (dnorm1*dnorm2)
 
 2873      if (cosphin >  m_one) cosphin =  m_one
 
 2874      if (cosphin < -m_one) cosphin = -m_one
 
 2876      phin = 
acos(cosphin)
 
 2879      point_1(1) = sfe(nsfe1)%x - sfe(ns)%x
 
 2880      point_1(2) = sfe(nsfe1)%y - sfe(ns)%y
 
 2881      point_1(3) = sfe(nsfe1)%z - sfe(ns)%z
 
 2883      dnorm1 = norm2(point_1)
 
 2885      if (abs(dnorm1) <= m_epsilon) dnorm1 = m_one
 
 2887      point_2(1) = pts(1,nn) - sfe(ns)%x
 
 2888      point_2(2) = pts(2,nn) - sfe(ns)%y
 
 2889      point_2(3) = pts(3,nn) - sfe(ns)%z
 
 2891      dnorm2 = norm2(point_2)
 
 2893      costn  = dot_product(point_1, point_2) / (dnorm1 * dnorm2)
 
 2894      sum1 = sum1 + phin * costn
 
 2905      if (nn > 1)   n0 = nn - 1
 
 2906      if (nn == 1)  n0 = nv
 
 2907      if (nn < nv)  n2 = nn + 1
 
 2908      if (nn == nv) n2 = 1
 
 2910      p1 = pts(:,n1) - ccc(:,n0)
 
 2911      p2 = pts(:,n0) - ccc(:,n0)
 
 2912      call vecp(p1, p2, p3, dnorm)
 
 2915      call vecp(p1, p2, p3, dnorm)
 
 2918      p1 = pts(:,n1) - ccc(:,n1)
 
 2919      p2 = pts(:,n2) - ccc(:,n1)
 
 2920      call vecp(p1, p2, p3, dnorm)
 
 2923      call vecp(p1, p2, p3, dnorm)
 
 2926      betan = 
acos(dot_product(u1, u2))
 
 2927      sum2 = sum2 + (m_pi - betan)
 
 2931    area = sfe(ns)%r*sfe(ns)%r*(tpi + sum1 - sum2)
 
 2937      pp(1) = pp(1) + (pts(1,ia) - sfe(ns)%x)
 
 2938      pp(2) = pp(2) + (pts(2,ia) - sfe(ns)%y)
 
 2939      pp(3) = pp(3) + (pts(3,ia) - sfe(ns)%z)
 
 2944    pp(1) = sfe(ns)%x + pp(1) * sfe(ns)%r / dnorm
 
 2945    pp(2) = sfe(ns)%y + pp(2) * sfe(ns)%r / dnorm
 
 2946    pp(3) = sfe(ns)%z + pp(3) * sfe(ns)%r / dnorm
 
 2949    pp1(1) = (pp(1) - sfe(ns)%x) / sfe(ns)%r
 
 2950    pp1(2) = (pp(2) - sfe(ns)%y) / sfe(ns)%r
 
 2951    pp1(3) = (pp(3) - sfe(ns)%z) / sfe(ns)%r
 
 2954    if (area < m_zero) area = m_zero
 
 2962  subroutine vecp(p1, p2, p3, dnorm)
 
 2963    real(real64), 
intent(in)  :: P1(:)
 
 2964    real(real64), 
intent(in)  :: P2(:)
 
 2965    real(real64), 
intent(out) :: P3(:)
 
 2966    real(real64), 
intent(out) :: dnorm
 
 2969    p3(1) = p1(2)*p2(3) - p1(3)*p2(2)
 
 2970    p3(2) = p1(3)*p2(1) - p1(1)*p2(3)
 
 2971    p3(3) = p1(1)*p2(2) - p1(2)*p2(1)
 
 2978    type(
pcm_t), 
intent(inout) :: pcm
 
 2980    integer :: asc_unit_test, asc_unit_test_sol, asc_unit_test_e, asc_unit_test_n, asc_unit_test_ext
 
 2985    if (pcm%solute .and. pcm%localf) 
then 
 2986      asc_unit_test     = io_open(pcm_dir//
'ASC.dat', pcm%namespace, action=
'write')
 
 2987      asc_unit_test_sol = io_open(pcm_dir//
'ASC_sol.dat', pcm%namespace, action=
'write')
 
 2988      asc_unit_test_e   = io_open(pcm_dir//
'ASC_e.dat', pcm%namespace, action=
'write')
 
 2989      asc_unit_test_n   = io_open(pcm_dir//
'ASC_n.dat', pcm%namespace, action=
'write')
 
 2990      asc_unit_test_ext = io_open(pcm_dir//
'ASC_ext.dat', pcm%namespace, action=
'write')
 
 2991      do ia = 1, pcm%n_tesserae
 
 2992        write(asc_unit_test,*)     pcm%tess(ia)%point, pcm%q_e(ia) + pcm%q_n(ia) + pcm%q_ext(ia), ia
 
 2993        write(asc_unit_test_sol,*) pcm%tess(ia)%point, pcm%q_e(ia) + pcm%q_n(ia), ia
 
 2994        write(asc_unit_test_e,*)   pcm%tess(ia)%point, pcm%q_e(ia), ia
 
 2995        write(asc_unit_test_n,*)   pcm%tess(ia)%point, pcm%q_n(ia), ia
 
 2996        write(asc_unit_test_ext,*) pcm%tess(ia)%point, pcm%q_ext(ia), ia
 
 2998      call io_close(asc_unit_test)
 
 2999      call io_close(asc_unit_test_sol)
 
 3000      call io_close(asc_unit_test_e)
 
 3001      call io_close(asc_unit_test_n)
 
 3002      call io_close(asc_unit_test_ext)
 
 3004    else if (pcm%solute .and. .not. pcm%localf) 
then 
 3005      asc_unit_test_sol = io_open(pcm_dir//
'ASC_sol.dat', pcm%namespace, action=
'write')
 
 3006      asc_unit_test_e   = io_open(pcm_dir//
'ASC_e.dat', pcm%namespace, action=
'write')
 
 3007      asc_unit_test_n   = io_open(pcm_dir//
'ASC_n.dat', pcm%namespace, action=
'write')
 
 3008      do ia = 1, pcm%n_tesserae
 
 3009        write(asc_unit_test_sol,*) pcm%tess(ia)%point, pcm%q_e(ia) + pcm%q_n(ia), ia
 
 3010        write(asc_unit_test_e,*)   pcm%tess(ia)%point, pcm%q_e(ia), ia
 
 3011        write(asc_unit_test_n,*)   pcm%tess(ia)%point, pcm%q_n(ia), ia
 
 3013      call io_close(asc_unit_test_sol)
 
 3014      call io_close(asc_unit_test_e)
 
 3015      call io_close(asc_unit_test_n)
 
 3017    else if (.not. pcm%solute .and. pcm%localf) 
then 
 3018      asc_unit_test_ext = io_open(pcm_dir//
'ASC_ext.dat', pcm%namespace, action=
'write')
 
 3019      do ia = 1, pcm%n_tesserae
 
 3020        write(asc_unit_test_ext,*) pcm%tess(ia)%point, pcm%q_ext(ia), ia
 
 3022      call io_close(asc_unit_test_ext)
 
 3025    safe_deallocate_a(pcm%spheres)
 
 3026    safe_deallocate_a(pcm%tess)
 
 3027    safe_deallocate_a(pcm%matrix)
 
 3028    if (.not. is_close(pcm%epsilon_infty, pcm%epsilon_0)) 
then 
 3029      safe_deallocate_a(pcm%matrix_d)
 
 3031    safe_deallocate_a(pcm%q_e)
 
 3032    safe_deallocate_a(pcm%q_e_in)
 
 3033    safe_deallocate_a(pcm%q_n)
 
 3034    safe_deallocate_a(pcm%v_e)
 
 3035    safe_deallocate_a(pcm%v_n)
 
 3036    safe_deallocate_a(pcm%v_e_rs)
 
 3037    safe_deallocate_a(pcm%v_n_rs)
 
 3038    if (pcm%localf) 
then 
 3039      safe_deallocate_a(pcm%matrix_lf)
 
 3040      if (.not. is_close(pcm%epsilon_infty, pcm%epsilon_0)) 
then 
 3041        safe_deallocate_a(pcm%matrix_lf_d)
 
 3044      safe_deallocate_a(pcm%q_ext)
 
 3045      safe_deallocate_a(pcm%q_ext_in)
 
 3046      safe_deallocate_a(pcm%v_ext)
 
 3047      safe_deallocate_a(pcm%v_ext_rs)
 
 3048      if (pcm%kick_is_present) 
then 
 3049        safe_deallocate_a(pcm%q_kick)
 
 3050        safe_deallocate_a(pcm%v_kick)
 
 3051        safe_deallocate_a(pcm%v_kick_rs)
 
 3056      safe_deallocate_a( pcm%rho_n)
 
 3057      safe_deallocate_a( pcm%rho_e)
 
 3058      if (pcm%localf) 
then 
 3059        safe_deallocate_a( pcm%rho_ext)
 
 3060        if (pcm%kick_is_present) 
then 
 3061          safe_deallocate_a( pcm%rho_kick)
 
 3066    if (pcm%tdlevel == 
pcm_td_eom) 
call pcm_eom_end()
 
 3068    if (mpi_grp_is_root(mpi_world)) 
call io_close(pcm%info_unit)
 
 3075  logical function pcm_update(this) 
result(update)
 
 3076    type(
pcm_t), 
intent(inout) :: this
 
 3078    this%iter = this%iter + 1
 
 3079    update = this%iter <= 6 .or. mod(this%iter, this%update_iter) == 0
 
 3085  real(real64) function 
pcm_get_vdw_radius(species, pcm_vdw_type, namespace)  result(vdw_r)
 
 3086    class(species_t),  
intent(in) :: species
 
 3087    integer,           
intent(in) :: pcm_vdw_type
 
 3088    type(namespace_t), 
intent(in) :: namespace
 
 3091    integer, 
parameter :: upto_xe = 54
 
 3092    real(real64)       :: vdw_radii(1:upto_xe)
 
 3095    data (vdw_radii(ia), ia=1, upto_xe)                                                                                  / &
 
 3097      1.001_real64,                                                                                            1.012_real64, &
 
 3099      0.825_real64, 1.408_real64,              1.485_real64, 2.000_real64, 1.583_real64, 1.500_real64, 1.287_real64, 1.243_real64, &
 
 3101      1.144_real64, 1.364_real64,              1.639_real64, 1.716_real64, 1.705_real64, 1.683_real64, 1.639_real64, 1.595_real64, &
 
 3103      1.485_real64, 1.474_real64,                                                                                            &
 
 3105      1.562_real64, 1.562_real64,                                                                    &
 
 3106      1.562_real64, 1.562_real64,                                                                    &
 
 3107      1.562_real64, 1.562_real64,                                                                    &
 
 3108      1.562_real64, 1.562_real64,                                                                    &
 
 3109      1.562_real64, 1.562_real64,                                                                  &
 
 3111      1.650_real64, 1.727_real64, 1.760_real64, 1.771_real64, 1.749_real64, 1.727_real64, &
 
 3113      1.628_real64, 1.606_real64, 1.639_real64, 1.639_real64,                                                                  &
 
 3114      1.639_real64, 1.639_real64,                                                            &
 
 3115      1.639_real64, 1.639_real64,                                                            &
 
 3116      1.639_real64, 1.639_real64,                                                                  &
 
 3117      1.639_real64, 1.639_real64,                                                                  &
 
 3119      2.672_real64, 1.804_real64, 1.881_real64, 1.892_real64, 1.892_real64, 1.881_real64  /
 
 3121    select case (pcm_vdw_type)
 
 3123      if (species%get_z() > upto_xe) 
then 
 3124        write(message(1),
'(a,a)') 
"The van der Waals radius is missing for element ", trim(species%get_label())
 
 3125        write(message(2),
'(a)') 
"Use PCMVdWRadii = pcm_vdw_species, for other vdw radii values" 
 3126        call messages_fatal(2, namespace=namespace)
 
 3128      ia = int(species%get_z())
 
 3129      vdw_r = vdw_radii(ia)*p_ang
 
 3132      vdw_r = species%get_vdw_radius()
 
 3133      if (vdw_r < m_zero) 
then 
 3134        call messages_write(
'The default vdW radius for species '//trim(species%get_label())//
':')
 
 3135        call messages_write(
' is not defined. ')
 
 3136        call messages_write(
' Add a positive vdW radius value in %Species block. ')
 
 3137        call messages_fatal(namespace=namespace)
 
 3146  subroutine pcm_dipole(mu_pcm, q_pcm, tess, n_tess)
 
 3147    real(real64),        
intent(out) :: mu_pcm(:)
 
 3148    real(real64),        
intent(in)  :: q_pcm(:)
 
 3149    integer,             
intent(in)  :: n_tess
 
 3150    type(pcm_tessera_t), 
intent(in)  :: tess(:)
 
 3158      mu_pcm = mu_pcm + q_pcm(ia) * tess(ia)%point
 
 3166  subroutine pcm_field(e_pcm, q_pcm, ref_point, tess, n_tess)
 
 3167    real(real64),        
intent(out) :: e_pcm(:)
 
 3168    real(real64),        
intent(in)  :: q_pcm(:)
 
 3169    integer,             
intent(in)  :: n_tess
 
 3170    type(pcm_tessera_t), 
intent(in)  :: tess(:)
 
 3171    real(real64),        
intent(in)  :: ref_point(1:3)
 
 3173    real(real64) :: diff(1:3)
 
 3174    real(real64) :: dist
 
 3182      diff = ref_point - tess(ia)%point
 
 3184      e_pcm = e_pcm + q_pcm(ia) * diff / dist**3
 
 3192  subroutine pcm_eps(pcm, eps, omega)
 
 3194    complex(real64), 
intent(out) :: eps
 
 3195    real(real64),    
intent(in)  :: omega
 
 3200      if (pcm%which_eps == pcm_debye_model) 
then 
 3202      else if (pcm%which_eps == pcm_drude_model) 
then 
 3218    type(namespace_t), 
intent(in)    :: namespace
 
 3223    call parse_variable(namespace, 
'PCMCalculation', .false., pcm%run_pcm)
 
 3224    call messages_print_with_emphasis(msg=
'PCM', namespace=namespace)
 
 3225    call parse_variable(namespace, 
'PCMLocalField', .false., pcm%localf)
 
 3226    call messages_print_var_value(
"PCMLocalField", pcm%localf, namespace=namespace)
 
 3227    if (pcm%localf) 
then 
 3228      call messages_experimental(
"PCM local field effects in the optical spectrum", namespace=namespace)
 
 3229      call messages_write(
'Beware of possible numerical errors in the optical spectrum due to PCM local field effects,')
 
 3230      call messages_new_line()
 
 3231      call messages_write(
'particularly, when static and high-frequency values of the dielectric functions are large')
 
 3232      call messages_write(
' (>~10 in units of the vacuum permittivity \epsilon_0).')
 
 3233      call messages_new_line()
 
 3234      call messages_write(
'However, PCM local field effects in the optical spectrum work well for polar or non-polar solvents')
 
 3235      call messages_new_line()
 
 3236      call messages_write(
'in the nonequilibrium or equation-of-motion TD-PCM propagation schemes.')
 
 3237      call messages_warning(namespace=namespace)
 
 3239    call parse_variable(namespace, 
'PCMTDLevel' , 
pcm_td_eq, pcm%tdlevel)
 
 3240    call messages_print_var_value(
"PCMTDLevel", pcm%tdlevel, namespace=namespace)
 
 3243    call parse_variable(namespace, 
'PCMStaticEpsilon' , m_one, pcm%deb%eps_0)
 
 3244    call messages_print_var_value(
"PCMStaticEpsilon", pcm%deb%eps_0, namespace=namespace)
 
 3246      call parse_variable(namespace, 
'PCMEpsilonModel', pcm_debye_model, pcm%which_eps)
 
 3247      call messages_print_var_value(
"PCMEpsilonModel", pcm%which_eps, namespace=namespace)
 
 3248      if (pcm%which_eps == pcm_debye_model) 
then 
 3249        call parse_variable(namespace, 
'PCMDynamicEpsilon', pcm%deb%eps_0, pcm%deb%eps_d)
 
 3250        call messages_print_var_value(
"PCMDynamicEpsilon", pcm%deb%eps_d, namespace=namespace)
 
 3251        call parse_variable(namespace, 
'PCMDebyeRelaxTime', m_zero, pcm%deb%tau)
 
 3252        call messages_print_var_value(
"PCMDebyeRelaxTime", pcm%deb%tau, namespace=namespace)
 
 3253      else if (pcm%which_eps == pcm_drude_model) 
then 
 3254        call parse_variable(namespace, 
'PCMDrudeLOmega', 
sqrt(m_one/(pcm%deb%eps_0-m_one)), pcm%drl%w0)
 
 3255        call messages_print_var_value(
"PCMDrudeLOmega", pcm%drl%w0, namespace=namespace)
 
 3256        call parse_variable(namespace, 
'PCMDrudeLDamping', m_zero, pcm%drl%gm)
 
 3257        call messages_print_var_value(
"PCMDrudeLDamping", pcm%drl%gm, namespace=namespace)
 
 3260      call parse_variable(namespace, 
'PCMDynamicEpsilon', pcm%deb%eps_0, pcm%deb%eps_d)
 
 3261      call messages_print_var_value(
"PCMDynamicEpsilon", pcm%deb%eps_d, namespace=namespace)
 
 3270    complex(real64),     
intent(out) :: eps
 
 3271    type(debye_param_t), 
intent(in) :: deb
 
 3272    real(real64),        
intent(in) :: omega
 
 3276    eps = deb%eps_d +  (deb%eps_0 - deb%eps_d)/(m_one + (omega*deb%tau)**2) + &
 
 3277      m_zi*omega*deb%tau*(deb%eps_0 - deb%eps_d)/(m_one + (omega*deb%tau)**2)
 
 3285    complex(real64),     
intent(out) :: eps
 
 3286    type(drude_param_t), 
intent(in)  :: drl
 
 3287    real(real64),        
intent(in)  :: omega
 
 3291    eps = m_one + (drl%w0**2 - omega**2)*drl%aa/((drl%w0**2 - omega**2)**2 + (omega*drl%gm)**2) + &
 
 3292      m_zi*omega*drl%gm*drl%aa/((drl%w0**2 - omega**2)**2 + (omega*drl%gm)**2)
 
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
 
double acos(double __x) __attribute__((__nothrow__
 
double exp(double __x) __attribute__((__nothrow__
 
double sin(double __x) __attribute__((__nothrow__
 
double sqrt(double __x) __attribute__((__nothrow__
 
double cos(double __x) __attribute__((__nothrow__
 
double floor(double __x) __attribute__((__nothrow__
 
type(debug_t), save, public debug
 
real(real64), parameter, public m_two
 
real(real64), parameter, public m_zero
 
real(real64), parameter, public m_four
 
real(real64), parameter, public p_a_b
some physical constants
 
real(real64), parameter, public m_pi
some mathematical constants
 
character(len= *), parameter, public pcm_dir
 
real(real64), parameter, public m_epsilon
 
real(real64), parameter, public m_half
 
real(real64), parameter, public m_one
 
This module implements the underlying real-space grid.
 
This module implements the index, used for the mesh points.
 
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
 
integer(int64) function, public io_function_fill_how(where)
Use this function to quickly plot functions for debugging purposes: call dio_function_output(io_funct...
 
subroutine, public io_close(iunit, grp)
 
subroutine, public io_mkdir(fname, namespace, parents)
 
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
 
This module is intended to contain "only mathematical" functions and procedures.
 
This module defines various routines, operating on mesh functions.
 
real(real64) function, public dmf_integrate(mesh, ff, mask, reduce)
Integrate a function on the mesh.
 
subroutine, public mesh_interpolation_init(this, mesh)
 
subroutine, public mesh_interpolation_end(this)
 
This module defines the meshes, which are used in Octopus.
 
integer function, public mesh_local_index_from_coords(mesh, ix)
This function returns the local index of the point for a given vector of integer coordinates.
 
integer(int64) function, public mesh_global_index_from_coords(mesh, ix)
This function returns the true global index of the point for a given vector of integer coordinates.
 
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
 
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
 
character(len=512), private msg
 
subroutine, public messages_warning(no_lines, all_nodes, namespace)
 
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
 
subroutine, public messages_new_line()
 
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
 
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
 
subroutine, public messages_input_error(namespace, var, details, row, column)
 
subroutine, public messages_experimental(name, namespace)
 
logical function mpi_grp_is_root(grp)
 
type(mpi_grp_t), public mpi_world
 
Some general things and nomenclature:
 
integer, parameter, public pcm_external_plus_kick
 
integer, parameter, public pcm_electrons
 
integer, parameter, public pcm_kick
 
integer, parameter, public pcm_nuclei
 
subroutine, public pcm_charges_propagation(q_t, pot_t, this_dt, this_cts_act, input_asc, this_eom, this_eps, namespace, this_deb, this_drl)
Driving subroutine for the Equation of Motion (EOM) propagation of the polarization charges within th...
 
integer, parameter, public pcm_external_potential
 
integer, parameter, public pcm_drude_model
 
subroutine, public pcm_eom_enough_initial(not_yet_called)
 
integer, parameter, public pcm_debye_model
 
subroutine, public pcm_calc_pot_rs(pcm, mesh, psolver, ions, v_h, v_ext, kick, time_present, kick_time)
 
subroutine pcm_eps_deb(eps, deb, omega)
 
logical function, public pcm_update(this)
Update pcm potential.
 
subroutine, public pcm_pot_rs(pcm, v_pcm, q_pcm, rho, mesh, psolver)
Generates the potential 'v_pcm' in real-space.
 
integer, parameter n_tess_sphere
minimum number of tesserae per sphere
 
subroutine pcm_eps_drl(eps, drl, omega)
 
logical function pcm_nn_in_mesh(pcm, mesh)
Check wether the nearest neighbor requested are in the mesh or not.
 
subroutine, public pcm_charge_density(pcm, q_pcm, q_pcm_tot, mesh, rho)
Generates the polarization charge density smearing the charge with a gaussian distribution on the mes...
 
real(real64) function, public pcm_get_vdw_radius(species, pcm_vdw_type, namespace)
get the vdw radius
 
integer, parameter pcm_dim_space
The resulting cavity is discretized by a set of tesserae defined in 3d.
 
subroutine pcm_matrix(eps, tess, n_tess, pcm_mat, localf)
Generates the PCM response matrix. J. Tomassi et al. Chem. Rev. 105, 2999 (2005).
 
subroutine, public pcm_dipole(mu_pcm, q_pcm, tess, n_tess)
Computes the dipole moment mu_pcm due to a distribution of charges q_pcm.
 
subroutine, public pcm_elect_energy(ions, pcm, E_int_ee, E_int_en, E_int_ne, E_int_nn, E_int_e_ext, E_int_n_ext)
Calculates the solute-solvent electrostatic interaction energy  .
 
real(real64) function s_mat_elem_i(tessi, tessj)
electrostatic Green function in vacuo:
 
subroutine d_i_matrix(n_tess, tess)
 
integer, parameter, public pcm_td_eom
 
subroutine, public pcm_charges(q_pcm, q_pcm_tot, v_cav, pcm_mat, n_tess, qtot_nominal, epsilon, renorm_charges, q_tot_tol, deltaQ)
Calculates the polarization charges at each tessera by using the response matrix 'pcm_mat',...
 
integer, parameter pcm_vdw_species
 
subroutine, public pcm_v_cav_li(v_cav, v_mesh, pcm, mesh)
Calculates the Hartree/external/kick potential at the tessera representative points by doing a 3D lin...
 
integer, parameter, public pcm_calc_poisson
 
subroutine gaubon(sfe, nv, ns, pts, ccc, pp, pp1, area, intsph)
Use the Gauss-Bonnet theorem to calculate the area of the tessera with vertices 'pts(3,...
 
real(real64), dimension(:,:), allocatable mat_gamess
PCM matrix formatted to be inputed to GAMESS.
 
subroutine pcm_poisson_sanity_check(pcm, mesh)
Check that all the required nearest neighbors are prensent in the mesh.
 
subroutine pcm_pot_rs_poisson(namespace, v_pcm, psolver, rho)
Generates the potential 'v_pcm' in real-space solving the poisson equation for rho.
 
real(real64) function d_mat_elem_i(tessi, tessj)
Gradient of the Green function in vacuo .
 
subroutine s_i_matrix(n_tess, tess)
 
subroutine pcm_pot_rs_direct(v_pcm, q_pcm, tess, n_tess, mesh, width_factor)
Generates the potential 'v_pcm' in real-space by direct sum.
 
logical gamess_benchmark
Decide to output pcm_matrix in a GAMESS format.
 
subroutine, public pcm_eps(pcm, eps, omega)
 
subroutine, public pcm_v_nuclei_cav(v_n_cav, ions, tess, n_tess)
Calculates the classical electrostatic potential geneated by the nuclei at the tesserae....
 
subroutine, public pcm_end(pcm)
 
integer, parameter pcm_vdw_optimized
 
subroutine, public pcm_field(e_pcm, q_pcm, ref_point, tess, n_tess)
Computes the field e_pcm at the reference point ref_point due to a distribution of charges q_pcm.
 
subroutine, public pcm_init(pcm, namespace, space, ions, grid, qtot, val_charge, external_potentials_present, kick_present)
Initializes the PCM calculation: reads the VdW molecular cavity and generates the PCM response matrix...
 
integer, parameter, public pcm_td_neq
 
subroutine vecp(p1, p2, p3, dnorm)
calculates the vectorial product p3 = p1 x p2
 
subroutine inter(sfe, p1, p2, p3, p4, ns, ia)
Finds the point 'p4', on the arc 'p1'-'p2' developed from 'p3', which is on the surface of sphere 'ns...
 
integer, parameter, public pcm_td_eq
 
subroutine subtessera(sfe, ns, nesf, nv, pts, ccc, pp, pp1, area)
find the uncovered region for each tessera and computes the area, the representative point (pp) and t...
 
subroutine, public pcm_min_input_parsing_for_spectrum(pcm, namespace)
 
subroutine cav_gen(tess_sphere, tess_min_distance, nesf, sfe, nts, cts, unit_pcminfo)
It builds the solute cavity surface and calculates the vertices, representative points and areas of t...
 
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel)
Calculates the Poisson equation. Given the density returns the corresponding potential.
 
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
 
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
 
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
 
This module defines the unit system, used for input and output.
 
type(unit_system_t), public units_out
 
type(unit_system_t), public units_inp
the units systems for reading and writing
 
type(unit_t), public unit_one
some special units required for particular quantities
 
Class implementing a box that is a union of spheres. We do this in a specific class instead of using ...
 
Description of the grid, containing information on derivatives, stencil, and symmetries.
 
Describes mesh distribution to nodes.
 
The cavity hosting the solute molecule is built from a set of interlocking spheres with optimized rad...