eplf/src/ezfio_interface.irp.f

165 lines
6.3 KiB
Fortran

BEGIN_SHELL [ /usr/bin/python ]
data = [ \
("nuclei_nucl_num" , "integer" , "" ),
("nuclei_nucl_charge" , "real" , "(nucl_num)" ),
("nuclei_nucl_coord" , "real" , "(nucl_num,3)" ),
("mo_basis_mo_coef" , "real" , "(ao_num,mo_tot_num)" ),
("mo_basis_mo_occ" , "real" , "(mo_tot_num)" ),
("electrons_elec_alpha_num" , "integer" , "" ),
("electrons_elec_beta_num" , "integer" , "" ),
("ao_basis_ao_num" , "integer" , "" ),
("ao_basis_ao_prim_num" , "integer" , "(ao_num)" ),
("ao_basis_ao_nucl" , "integer" , "(ao_num)" ),
("ao_basis_ao_power" , "integer" , "(ao_num,3)" ),
("ao_basis_ao_expo" , "real" , "(ao_num,ao_prim_num_max)" ),
("ao_basis_ao_coef" , "real" , "(ao_num,ao_prim_num_max)" ),
("determinants_det_num" , "integer" , "" ),
("determinants_det_coef" , "real" , "(det_num)" ),
("determinants_det_occ" , "integer" , "(elec_alpha_num-mo_closed_num,det_num,2)" ),
("grid_point_num" , "integer" , "(3)" ),
("grid_step_size" , "real" , "(3)" ),
("grid_origin" , "real" , "(3)" ),
("grid_opposite" , "real" , "(3)" ),
("compute_eplf" , "logical" , "" ),
("compute_eplf_grad" , "logical" , "" ),
("compute_eplf_lapl" , "logical" , "" ),
("compute_elf" , "logical" , "" ),
("compute_elf_grad" , "logical" , "" ),
("compute_elf_lapl" , "logical" , "" ),
("compute_density" , "logical" , "" ),
("compute_density_grad" , "logical" , "" ),
("compute_density_lapl" , "logical" , "" ),
("compute_elf_partition" , "logical" , "" ),
("compute_density_partition" , "logical" , "" ),
("compute_eplf_partition" , "logical" , "" ),
("grid_data_eplf" , "real" , "(grid_x_num,grid_y_num,grid_z_num)" ),
("grid_data_eplf_grad" , "real" , "(grid_x_num,grid_y_num,grid_z_num,4)" ),
("grid_data_eplf_lapl" , "real" , "(grid_x_num,grid_y_num,grid_z_num)" ),
("grid_data_elf" , "real" , "(grid_x_num,grid_y_num,grid_z_num)" ),
("grid_data_elf_grad" , "real" , "(grid_x_num,grid_y_num,grid_z_num,4)" ),
("grid_data_elf_lapl" , "real" , "(grid_x_num,grid_y_num,grid_z_num)" ),
("grid_data_density" , "real" , "(grid_x_num,grid_y_num,grid_z_num)" ),
("grid_data_density_grad" , "real" , "(grid_x_num,grid_y_num,grid_z_num,4)" ),
("grid_data_density_lapl" , "real" , "(grid_x_num,grid_y_num,grid_z_num)" ),
("grid_data_elf_partition" , "real" , "(grid_x_num,grid_y_num,grid_z_num)" ),
("grid_data_eplf_partition" , "real" , "(grid_x_num,grid_y_num,grid_z_num)" ),
("grid_data_density_partition" , "real" , "(grid_x_num,grid_y_num,grid_z_num)" ),
("density_matrix_one" , "real" , "(mo_active_num,mo_active_num,2)" ),
("density_matrix_two_num" , "integer" , "" ),
("density_matrix_two_indice" , "integer" , "(4,two_e_density_num_max)" ),
("density_matrix_two_value" , "real" , "(2,two_e_density_num_max)" ),
]
data_no_set = [\
("mo_basis_mo_tot_num" , "integer" , "" ),
("mo_basis_mo_active_num" , "integer" , "" ),
("mo_basis_mo_closed_num" , "integer" , "" ),
]
def do_subst(t0,d):
t = t0
t = t.replace("$X",d[0])
t = t.replace("$T",d[1])
t = t.replace("$D",d[2])
if d[1].startswith("character"):
size = d[1].split("*")[1][1:-1]
u = "character"
elif d[1].startswith("double precision"):
u = d[1].replace(" ","_")
size = "1"
elif "*" in d[1]:
size = "1"
u = d[1].replace("*","")
else:
size = "1"
u = d[1]
t = t.replace("$U",u)
if d[2] == "":
t = t.replace("$S",size)
else:
if size == "1":
t = t.replace("$S","size(res)")
else:
t = t.replace("$S","%s*size(res)"%(size))
print t
t0 = """
subroutine get_$X(res)
implicit none
IRP_IF MPI
include 'mpif.h'
IRP_ENDIF
$T :: res$D
integer :: ierr
logical :: exists
!$OMP CRITICAL (ezfio_critical)
PROVIDE ezfio_filename
if (mpi_master) then
call ezfio_has_$X(exists)
if (exists) then
call ezfio_get_$X(res)
else
call ezfio_set_$X(res)
endif
endif
IRP_IF MPI
call MPI_BCAST(res,$S,MPI_$U,0,MPI_COMM_WORLD,ierr)
if (ierr /= MPI_SUCCESS) then
call abrt(irp_here,'Unable to broadcast $X')
endif
IRP_ENDIF
!$OMP END CRITICAL (ezfio_critical)
end
subroutine set_$X(res)
implicit none
$T :: res$D
integer :: ierr
logical :: exists
PROVIDE ezfio_filename
!$OMP CRITICAL (ezfio_critical)
if (mpi_master) then
call ezfio_set_$X(res)
endif
!$OMP END CRITICAL (ezfio_critical)
end
"""
t1 = """
subroutine get_$X(res)
implicit none
IRP_IF MPI
include 'mpif.h'
IRP_ENDIF
$T :: res$D
integer :: ierr
!$OMP CRITICAL (ezfio_critical)
PROVIDE ezfio_filename
if (mpi_master) then
call ezfio_get_$X(res)
endif
IRP_IF MPI
call MPI_BCAST(res,$S,MPI_$U,0,MPI_COMM_WORLD,ierr)
if (ierr /= MPI_SUCCESS) then
call abrt(irp_here,'Unable to broadcast $X')
endif
IRP_ENDIF
!$OMP END CRITICAL (ezfio_critical)
end
"""
for d in data:
do_subst(t0,d)
for d in data_no_set:
do_subst(t1,d)
END_SHELL