1
0
mirror of https://github.com/TREX-CoE/trexio.git synced 2024-11-03 20:54:07 +01:00

Merge pull request #89 from TREX-CoE/walkers

Electron coordinates
This commit is contained in:
Evgeny Posenitskiy 2022-04-13 16:38:58 +02:00 committed by GitHub
commit dc783bc1e0
No known key found for this signature in database
GPG Key ID: 4AEE18F83AFDEB23
7 changed files with 266 additions and 23 deletions

View File

@ -117,7 +117,8 @@ TESTS_C = \
tests/io_dset_str_text \
tests/delete_group_text \
tests/overwrite_all_text \
tests/io_all
tests/io_all \
tests/pre_close
if HAVE_HDF5
TESTS_C += \
@ -252,10 +253,15 @@ $(pytrexio_py): $(pytrexio_c)
# Build Python module and C wrapper code for TREXIO using SWIG
# [?] swig -python -threads pytrexio.i ----> Add thread support for all the interface
$(pytrexio_c): $(ORG_FILES) $(GENERATOR_FILES) $(trexio_h) $(pytrexio_i) $(numpy_i)
cp $(trexio_h) src/
@if [[ $(SWIG).x != ".x" ]] ; then \
cp $(trexio_h) src/ ; \
cd src/ && \
$(SWIG) -python -py3 -o pytrexio_wrap.c pytrexio.i
$(RM) -- src/trexio.h
$(SWIG) -python -py3 -o pytrexio_wrap.c pytrexio.i ; \
$(RM) -- src/trexio.h ;\
else echo "Error: SWIG is not installed" ; \
exit 1 ; \
fi
$(numpy_i):
wget https://raw.githubusercontent.com/numpy/numpy/main/tools/swig/numpy.i -O $(numpy_i)

View File

@ -8,7 +8,7 @@
#+INFOJS_OPT: toc:t mouse:underline path:org-info.js
#+HTML_HEAD: <link rel="stylesheet" title="Standard" href="trexio.css" type="text/css" />
#+STARTUP: align fold nodlcheck hidestars oddeven lognotestate
#+STARTUP: align nodlcheck hidestars oddeven lognotestate
#+AUTHOR: TREX-CoE
#+LANGUAGE: en

View File

@ -864,19 +864,17 @@ trexio_open(const char* file_name, const char mode,
/* Try to determine the applicable backend if the back_end argument is TREXIO_AUTO */
if (back_end == TREXIO_AUTO && mode == 'r') {
#ifdef HAVE_HDF5
trexio_exit_code rc_text, rc_hdf5;
/* Check if the TREXIO file exists and if it is a directory */
rc_text = trexio_text_inquire(file_name);
trexio_exit_code rc_text = trexio_text_inquire(file_name);
if (rc_text == TREXIO_SUCCESS) {
back_end_local = TREXIO_TEXT;
} else {
/* If not, check if it is an HDF5 file */
rc_hdf5 = trexio_hdf5_inquire(file_name);
trexio_exit_code rc_hdf5 = trexio_hdf5_inquire(file_name);
if (rc_hdf5 == TREXIO_SUCCESS) {
back_end_local = TREXIO_HDF5;
} else {
/* File is neither a directory nor an HDF5 file -> return an error */
back_end_local = -1;
if (rc_open != NULL) *rc_open = TREXIO_FILE_ERROR;
return NULL;
}
@ -1211,6 +1209,12 @@ trexio_close (trexio_t* file)
assert(file->back_end < TREXIO_INVALID_BACK_END);
/* Things to be done before the closing the file in the back-end */
rc = trexio_pre_close(file);
if (rc != TREXIO_SUCCESS) {
return rc;
}
/* Terminate the back end */
switch (file->back_end) {
@ -1388,6 +1392,82 @@ def _inquire(file_name: str) -> bool:
raise Error(rc)
#+end_src
** Tasks to be done before closing
#+begin_src c :tangle trexio_private.h :exports none
trexio_exit_code trexio_pre_close(trexio_t* file);
#+end_src
#+begin_src c :tangle prefix_front.c
trexio_exit_code
trexio_pre_close (trexio_t* file)
{
if (file == NULL) return TREXIO_FILE_ERROR;
{ /* Up-spin and down-spin electrons */
trexio_exit_code rc;
int32_t nup, ndn, nelec;
bool has_up = (trexio_has_electron_up_num(file) == TREXIO_SUCCESS);
bool has_dn = (trexio_has_electron_dn_num(file) == TREXIO_SUCCESS);
bool has_updn = (trexio_has_electron_num(file) == TREXIO_SUCCESS);
if (has_updn && has_up && has_dn) {
rc = trexio_read_electron_up_num(file, &nup);
if (rc != TREXIO_SUCCESS) return rc;
rc = trexio_read_electron_dn_num(file, &ndn);
if (rc != TREXIO_SUCCESS) return rc;
rc = trexio_read_electron_num(file, &nelec);
if (rc != TREXIO_SUCCESS) return rc;
if (nelec != nup + ndn) {
nelec = nup + ndn;
rc = trexio_write_electron_num(file, nelec);
if (rc != TREXIO_SUCCESS) return rc;
}
} else if (has_up && has_dn) {
rc = trexio_read_electron_up_num(file, &nup);
if (rc != TREXIO_SUCCESS) return rc;
rc = trexio_read_electron_dn_num(file, &ndn);
if (rc != TREXIO_SUCCESS) return rc;
nelec = nup + ndn;
rc = trexio_write_electron_num(file, nelec);
if (rc != TREXIO_SUCCESS) return rc;
} else if (has_up) {
rc = trexio_read_electron_up_num(file, &nup);
if (rc != TREXIO_SUCCESS) return rc;
ndn = 0;
rc = trexio_write_electron_dn_num(file, ndn);
if (rc != TREXIO_SUCCESS) return rc;
nelec = nup;
rc = trexio_write_electron_num(file, nelec);
if (rc != TREXIO_SUCCESS) return rc;
} else if (has_dn) {
rc = trexio_read_electron_dn_num(file, &ndn);
if (rc != TREXIO_SUCCESS) return rc;
nup = 0;
rc = trexio_write_electron_up_num(file, nup);
if (rc != TREXIO_SUCCESS) return rc;
nelec = ndn;
rc = trexio_write_electron_num(file, nelec);
if (rc != TREXIO_SUCCESS) return rc;
}
}
return TREXIO_SUCCESS;
}
#+end_src
* Templates for front end
** Description

View File

@ -439,7 +439,7 @@ trexio_hdf5_write_$group_dset$ (trexio_t* const file,
trexio_hdf5_t* f = (trexio_hdf5_t*) file;
hid_t index_dtype;
void* index_p = NULL;
const void* index_p;
uint64_t size_ranked = (uint64_t) size * $group_dset_rank$;
/* Determine the optimal type for storing indices depending on the size_max (usually mo_num or ao_num) */
if (size_max < UINT8_MAX) {
@ -459,7 +459,7 @@ trexio_hdf5_write_$group_dset$ (trexio_t* const file,
index_p = index;
index_dtype = H5T_NATIVE_UINT16;
} else {
index_p = (int32_t*) index_sparse;
index_p = (const int32_t*) index_sparse;
index_dtype = H5T_NATIVE_INT32;
}

View File

@ -13,6 +13,7 @@ set(Tests_text
io_str_text
delete_group_text
overwrite_all_text
pre_close
)
if(ENABLE_HDF5)

124
tests/pre_close.c Normal file
View File

@ -0,0 +1,124 @@
#include "trexio.h"
#include <assert.h>
#include <stdio.h>
#include <stdlib.h>
#include <stdint.h>
#define TEST_BACKEND TREXIO_TEXT
#define TREXIO_FILE "test_dset_sparse.dir"
#define RM_COMMAND "rm -rf " TREXIO_FILE
static int test_pre_close_1 (const char* file_name, const back_end_t backend)
{
/* Check if nelec = nup + ndn */
trexio_t* file = NULL;
trexio_exit_code rc;
/*================= START OF TEST ==================*/
// open file in 'write' mode
file = trexio_open(file_name, 'w', backend, &rc);
assert (file != NULL);
assert (rc == TREXIO_SUCCESS);
// write parameters
int32_t nup = 4;
int32_t ndn = 3;
int32_t nelec = 0;
rc = trexio_write_electron_up_num(file, nup);
assert (rc == TREXIO_SUCCESS);
rc = trexio_write_electron_dn_num(file, ndn);
assert (rc == TREXIO_SUCCESS);
// close file
rc = trexio_close(file);
assert (rc == TREXIO_SUCCESS);
// re-open file
file = trexio_open(file_name, 'r', backend, &rc);
assert (file != NULL);
assert (rc == TREXIO_SUCCESS);
rc = trexio_read_electron_num(file, &nelec);
assert (rc == TREXIO_SUCCESS);
printf("nup : %d\n", nup);
printf("ndn : %d\n", ndn);
printf("nelec: %d\n", nelec);
assert (nelec == nup + ndn);
// close file
rc = trexio_close(file);
assert (rc == TREXIO_SUCCESS);
/*================= END OF TEST ==================*/
return 0;
}
static int test_pre_close_2 (const char* file_name, const back_end_t backend)
{
/* Check if nelec = nup */
trexio_t* file = NULL;
trexio_exit_code rc;
/*================= START OF TEST ==================*/
// open file in 'write' mode
file = trexio_open(file_name, 'w', backend, &rc);
assert (file != NULL);
assert (rc == TREXIO_SUCCESS);
// write parameters
int32_t nup = 4;
int32_t nelec = 0;
rc = trexio_write_electron_up_num(file, nup);
assert (rc == TREXIO_SUCCESS);
// close file
rc = trexio_close(file);
assert (rc == TREXIO_SUCCESS);
// re-open file
file = trexio_open(file_name, 'r', backend, &rc);
assert (file != NULL);
assert (rc == TREXIO_SUCCESS);
rc = trexio_read_electron_num(file, &nelec);
assert (rc == TREXIO_SUCCESS);
assert (nelec == nup);
// close file
rc = trexio_close(file);
assert (rc == TREXIO_SUCCESS);
/*================= END OF TEST ==================*/
return 0;
}
int main()
{
/*============== Test launcher ================*/
int rc;
rc = system(RM_COMMAND);
assert (rc == 0);
test_pre_close_1 (TREXIO_FILE, TEST_BACKEND);
rc = system(RM_COMMAND);
test_pre_close_2 (TREXIO_FILE, TEST_BACKEND);
rc = system(RM_COMMAND);
assert (rc == 0);
return 0;
}

View File

@ -5,14 +5,10 @@
This page contains information about the general structure of the
TREXIO library. The source code of the library can be automatically
generated based on the contents of the ~trex.json~ configuration file,
which itself is compiled from different sections (groups) presented
which itself is generated from different sections (groups) presented
below.
For more information about the automatic generation on the source code
or regarding possible modifications, please contact the TREXIO
developers.
All quantities are saved in TREXIO file in atomic units. The
All quantities are saved in TREXIO files in atomic units. The
dimensions of the arrays in the tables below are given in column-major
order (as in Fortran), and the ordering of the dimensions is reversed
in the produced ~trex.json~ configuration file as the library is
@ -21,19 +17,21 @@ written in C.
TREXIO currently supports ~int~, ~float~ and ~str~ types for both
single attributes and arrays. Note, that some attributes might have
~dim~ type (e.g. ~num~ of the ~nucleus~ group). This type is treated
exactly the same as ~int~ with the only difference that ~dim~
exactly in the same way as ~int~ with the only difference that ~dim~
variables cannot be negative. This additional constraint is required
because ~dim~ attributes are used internally to allocate memory and to
check array boundaries in the memory-safe API. Most of the times, the
~dim~ variables contain the ~num~ suffix.
In Fortran, the arrays are 1-based and in most other languages the
arrays are 0-based. Hence, we introduce the ~index~ type which is an
In Fortran, arrays are 1-based and in most other languages the
arrays are 0-based. Hence, we introduce the ~index~ type which is a
1-based ~int~ in the Fortran interface and 0-based otherwise.
For sparse data structures such as electron replusion integrals,
the data can be too large to fit in memory and the data needs to be
fetched using multiple function calls to perform I/O on buffers.
For more information on how to read/write sparse data structures, see
the [[./examples.html][examples]].
#+begin_src python :tangle trex.json :exports none
@ -88,6 +86,7 @@ fetched using multiple function calls to perform I/O on buffers.
#+NAME:electron
| Variable | Type | Dimensions | Description |
|----------+-------+------------+-------------------------------------|
| ~num~ | ~dim~ | | Number of electrons |
| ~up_num~ | ~int~ | | Number of \uparrow-spin electrons |
| ~dn_num~ | ~int~ | | Number of \downarrow-spin electrons |
@ -96,7 +95,8 @@ fetched using multiple function calls to perform I/O on buffers.
:results:
#+begin_src python :tangle trex.json
"electron": {
"up_num" : [ "int", [] ]
"num" : [ "dim", [] ]
, "up_num" : [ "int", [] ]
, "dn_num" : [ "int", [] ]
} ,
#+end_src
@ -689,7 +689,7 @@ prim_factor =
| ~2e_updn~ | ~float sparse~ | ~(mo.num, mo.num, mo.num, mo.num)~ | \uparrow\downarrow component of the two-body reduced density matrix |
| ~2e_dnup~ | ~float sparse~ | ~(mo.num, mo.num, mo.num, mo.num)~ | \downarrow\uparrow component of the two-body reduced density matrix |
#+CALL: json(data=rdm, title="rdm", last=1)
#+CALL: json(data=rdm, title="rdm")
#+RESULTS:
:results:
@ -703,10 +703,42 @@ prim_factor =
, "2e_dndn" : [ "float sparse", [ "mo.num", "mo.num", "mo.num", "mo.num" ] ]
, "2e_updn" : [ "float sparse", [ "mo.num", "mo.num", "mo.num", "mo.num" ] ]
, "2e_dnup" : [ "float sparse", [ "mo.num", "mo.num", "mo.num", "mo.num" ] ]
}
} ,
#+end_src
:end:
* Quantum Monte Carlo data (qmc group)
In quantum Monte Carlo calculations, the wave function is evaluated
at points of the 3N-dimensional space. Some algorithms require multiple
independent /walkers/, so it is possible to store multiple coordinates,
as well as some quantities evaluated at those points.
By convention, the electron coordinates contain first all the electrons
of $\uparrow$-spin and then all the $\downarrow$-spin.
#+name: qmc
| Variable | Type | Dimensions | Description |
|----------+---------+------------------------------+---------------------------------------|
| ~num~ | ~dim~ | | Number of 3N-dimensional points |
| ~point~ | ~float~ | ~(3, electron.num, qmc.num)~ | 3N-dimensional points |
| ~psi~ | ~float~ | ~(qmc.num)~ | Wave function evaluated at the points |
| ~e_loc~ | ~float~ | ~(qmc.num)~ | Local energy evaluated at the points |
#+CALL: json(data=qmc, title="qmc", last=1)
#+RESULTS:
:results:
#+begin_src python :tangle trex.json
"qmc": {
"num" : [ "dim" , [] ]
, "point" : [ "float", [ "qmc.num", "electron.num", "3" ] ]
, "psi" : [ "float", [ "qmc.num" ] ]
, "e_loc" : [ "float", [ "qmc.num" ] ]
}
#+end_src
:end:
* Appendix :noexport:
** Python script from table to json