1 /*****************************************************************************
2  *  CP2K: A general program to perform molecular dynamics simulations        *
3  *  Copyright (C) 2000 - 2020  CP2K developers group                         *
4  *****************************************************************************/
5 
6 #include <stdbool.h>
7 
8 /**
9  * \brief Definitions for the functions exported in libcp2k.F
10  * \author Mohammad Hossein Bani-Hashemian
11  */
12 
13 #ifndef LIBCP2K_H
14 #define LIBCP2K_H
15 
16 #ifdef __cplusplus
17 extern "C" {
18 #endif
19 
20 typedef int force_env_t;
21 
22 /** \brief Get the CP2K version string
23  * \param version_str The buffer to write the version string into
24  * \param str_length The size of the buffer (must be large enough)
25  */
26 void cp2k_get_version(char* version_str, int str_length);
27 
28 /** \brief Initialize CP2K and MPI
29  * \warning You are supposed to call cp2k_finalize() before terminating the program.
30  */
31 void cp2k_init();
32 
33 /** \brief Initialize CP2K without initializing MPI
34  * \warning You are supposed to call cp2k_finalize() before terminating the program.
35  */
36 void cp2k_init_without_mpi();
37 
38 /** \brief Finalize CP2K and MPI
39  */
40 void cp2k_finalize();
41 
42 /** \brief Finalize CP2K and without finalizing MPI
43  */
44 void cp2k_finalize_without_mpi();
45 
46 /** \brief Create a new force environment
47  * \param new_force_env the created force environment
48  * \param input_file_path Path to a CP2K input file
49  * \param output_file_path Path to a file where CP2K is going to append its output (created if non-existent)
50  * \warning You are supposed to call cp2k_destroy_force_env() to cleanup, before cp2k_finalize().
51  */
52 void cp2k_create_force_env(force_env_t* new_force_env, const char* input_file_path, const char* output_file_path);
53 
54 /** \brief Create a new force environment (custom managed MPI)
55  * \param new_force_env the created force environment
56  * \param input_file_path Path to a CP2K input file
57  * \param output_file_path Path to a file where CP2K is going to write its output.
58  *                         Will be created if not existent, otherwise appended.
59  * \param mpi_comm Fortran MPI communicator if MPI is not managed by CP2K
60  * \warning You are supposed to call cp2k_destroy_force_env() to cleanup, before cp2k_finalize().
61  */
62 void cp2k_create_force_env_comm(force_env_t* new_force_env, const char* input_file_path, const char* output_file_path, int mpi_comm);
63 
64 /** \brief Destroy/cleanup a force environment
65  * \param force_env the force environment
66  */
67 void cp2k_destroy_force_env(force_env_t force_env);
68 
69 /** \brief Set positions of the particles
70  * \param force_env the force environment
71  * \param new_pos Array containing the new positions of the particles
72  * \param n_el Size of the new_pos array
73  */
74 void cp2k_set_positions(force_env_t force_env, const double* new_pos, int n_el);
75 
76 /** \brief Set velocity of the particles
77  * \param force_env the force environment
78  * \param new_vel Array containing the new velocities of the particles
79  * \param n_el Size of the new_vel array
80  */
81 void cp2k_set_velocities(force_env_t force_env, const double* new_vel, int n_el);
82 
83 /** \brief Get an arbitrary result as 1D array from CP2K
84  * \param force_env the force environment
85  * \param description The string tag of the result
86  * \param results Pre-allocated array
87  * \param n_el size of the results array
88  */
89 void cp2k_get_result(force_env_t force_env, const char* description, double* result, int n_el);
90 
91 /** \brief Get the number of atoms
92  * \param force_env the force environment
93  * \param natom The number of atoms
94  */
95 void cp2k_get_natom(force_env_t force_env, int* natom);
96 
97 /** \brief Get the number of particles
98  * \param force_env the force environment
99  * \param nparticle The number of particles
100  */
101 void cp2k_get_nparticle(force_env_t force_env, int* nparticle);
102 
103 /** \brief Get the positions of the particles
104  * \param force_env the force environment
105  * \param pos Pre-allocated array of at least 3*nparticle elements. Use cp2k_get_nparticle() to get the number of particles.
106  * \param n_el Size of the force array
107  */
108 void cp2k_get_positions(force_env_t force_env, double* pos, int n_el);
109 
110 /** \brief Get the forces for the particles
111  * \param force_env the force environment
112  * \param force Pre-allocated array of at least 3*nparticle elements. Use cp2k_get_nparticle() to get the number of particles.
113  * \param n_el Size of the force array
114  */
115 void cp2k_get_forces(force_env_t force_env, double* force, int n_el);
116 
117 /** \brief Get the potential energy of the system
118  * \param force_env the force environment
119  * \param e_pot The potential energy
120  */
121 void cp2k_get_potential_energy(force_env_t force_env, double* e_pot);
122 
123 /** \brief Calculate energy and forces of the system
124  * \param force_env the force environment
125  */
126 void cp2k_calc_energy_force(force_env_t force_env);
127 
128 /** \brief Calculate only the energy of the system
129  * \param force_env the force environment
130  */
131 void cp2k_calc_energy(force_env_t force_env);
132 
133 /** \brief Make a CP2K run with the given input file
134  * \param input_file_path Path to a CP2K input file
135  * \param output_file_path Path to a file where CP2K is going to append its output (created if non-existent)
136  */
137 void cp2k_run_input(const char* input_file_path, const char* output_file_path);
138 
139 /** \brief Make a CP2K run with the given input file (custom managed MPI)
140  * \param input_file_path Path to a CP2K input file
141  * \param output_file_path Path to a file where CP2K is going to append its output (created if non-existent)
142  * \param mpi_comm Fortran MPI communicator if MPI is not managed by CP2K
143  */
144 void cp2k_run_input_comm(const char* input_file_path, const char* output_file_path, int mpi_comm);
145 
146 /** \brief Transport parameters read from a CP2K input file
147  *
148  * This definition matches the respective type definition in the transport_env_types module
149  */
150 typedef struct {
151     int     n_occ;
152     int     n_atoms;
153     double  energy_diff;
154     double  evoltfactor;
155     double  e_charge;
156     double  boltzmann;
157     double  h_bar;
158     int     iscf;
159     int     method;
160     int     qt_formalism;
161     int     injection_method;
162     int     rlaxis_integration_method;
163     int     linear_solver;
164     int     matrixinv_method;
165     int     transport_neutral;
166     int     num_pole;
167     int     ordering;
168     int     row_ordering;
169     int     verbosity;
170     int     pexsi_np_symb_fact;
171     int     n_kpoint;
172     int     num_interval;
173     int     num_contacts;
174     int     stride_contacts;
175     int     tasks_per_energy_point;
176     int     tasks_per_pole;
177     int     gpus_per_point;
178     int     n_points_beyn;
179     int     ncrc_beyn;
180     int     tasks_per_integration_point;
181     int     n_points_inv;
182     int     cutout[2];
183     double  colzero_threshold;
184     double  eps_limit;
185     double  eps_limit_cc;
186     double  eps_decay;
187     double  eps_singularity_curvatures;
188     double  eps_mu;
189     double  eps_eigval_degen;
190     double  eps_fermi;
191     double  energy_interval;
192     double  min_interval;
193     double  temperature;
194     double  dens_mixing;
195     double  n_rand_beyn;
196     double  n_rand_cc_beyn;
197     double  svd_cutoff;
198     int*    contacts_data;
199     int*    nsgf;
200     double* zeff;
201     bool    obc_equilibrium;
202     bool    extra_scf;
203 } cp2k_transport_parameters;
204 
205 /** \brief CP2K's C-interoperable CSR matrix
206  *
207  * This definition matches the respective type definition in the transport_env_types module
208  */
209 typedef struct {
210     int     nrows_total;
211     int     ncols_total;
212     int     nze_total;
213     int     nze_local;
214     int     nrows_local;
215     int     data_type;
216     int     first_row;
217     int*    rowptr_local;
218     int*    colind_local;
219     int*    nzerow_local;
220     double* nzvals_local;
221 } cp2k_csr_interop_type;
222 
223 /** \brief Function pointer type for the externally evaluated density matrix
224  *
225  * Function pointer type pointing to a C routine that takes the S and H matrices as input and outputs a P matrix.
226  *
227  * Function definition example:
228  * \code{.c}
229  * void c_scf_method(
230  *     cp2k_transport_parameters cp2k_transport_params,
231  *     cp2k_csr_interop_type S,
232  *     cp2k_csr_interop_type KS,
233  *     cp2k_csr_interop_type* P,
234  *     cp2k_csr_interop_type* PImag
235  *     );
236  * \endcode
237  * \sa cp2k_transport_parameters, cp2k_csr_interop_type
238  */
239 typedef void (*ext_method_callback_f_ptr) (
240     cp2k_transport_parameters, // Transport parameters
241     cp2k_csr_interop_type,  // S-Matrix
242     cp2k_csr_interop_type,  // H-Matrix
243     cp2k_csr_interop_type*, // P-Matrix
244     cp2k_csr_interop_type*  // PImag-Matrix
245     );
246 
247 /** \brief Set the function callback for the externally evaluated density matrix
248  */
249 void cp2k_transport_set_callback(force_env_t force_env, ext_method_callback_f_ptr func);
250 
251 /** \brief Get the number of molecular orbitals in the active space
252  * \param force_env the force environment
253  * \returns The number of elements or -1 if unavailable
254  */
255 int cp2k_active_space_get_mo_count(force_env_t force_env);
256 
257 /** \brief Get the Fock submatrix for the active space
258  * \param force_env the force environment
259  * \param buf Pre-allocated array of at least mo_count^2 elements.
260  *            Use `cp2k_active_space_get_mo_count()` to get the number of molecular orbitals.
261  * \param buf_len Size of the buf array
262  * \returns The number of elements written to buf or -1 if unavailable
263  */
264 long int cp2k_active_space_get_fock_sub(force_env_t force_env, double* buf, long int buf_len);
265 
266 /** \brief Get the number of non-zero elements in the ERI matrix
267  * \param force_env the force environment
268  * \returns The number of elements or -1 if unavailable
269  */
270 long int cp2k_active_space_get_eri_nze_count(force_env_t force_env);
271 
272 /** \brief Get the non-zero elements of the ERI matrix
273  *
274  * The buf_coords will contain the coordinates in the format `[i1, j1, k1, l1, i2, j2, k2, l2, ... ]`.
275  * \param force_env the force environment
276  * \param buf_coords Pre-allocated array of at least 4*nze_count elements.
277  *                   Use `cp2k_active_space_get_eri_nze_count()` to get the number of non-zero elements.
278  * \param buf_coords_len Size of the buf_coords array
279  * \param buf_values Pre-allocated array of at least nze_count elements.
280  *                   Use `cp2k_active_space_get_eri_nze_count()` to get the number of non-zero elements.
281  * \param buf_values_len Size of the buf_values array
282  * \returns The number of elements written to buf_values or -1 if unavailable
283  */
284 int cp2k_active_space_get_eri(force_env_t force_env,
285                               int* buf_coords, long int buf_coords_len,
286                               double* buf_values, long int buf_values_len);
287 
288 #ifdef __cplusplus
289 }
290 #endif
291 
292 #endif
293 /* vim: set ts=4 sw=4 tw=0 :*/
294