1
0
mirror of https://github.com/TREX-CoE/qmckl.git synced 2024-07-22 18:57:40 +02:00

1-body jastrow dimensioned by nucl type

This commit is contained in:
Anthony Scemama 2022-11-17 16:54:49 +01:00
parent 5a6d064ce0
commit da2d8f6250
3 changed files with 62 additions and 50 deletions

View File

@ -2658,11 +2658,10 @@ qmckl_exit_code qmckl_finalize_basis(qmckl_context context) {
#ifdef HAVE_HPC #ifdef HAVE_HPC
rc = qmckl_finalize_basis_hpc(context); rc = qmckl_finalize_basis_hpc(context);
#else if (rc != QMCKL_SUCCESS) return rc;
rc = QMCKL_SUCCESS;
#endif #endif
return rc; return qmckl_context_touch(context);
} }
#+end_src #+end_src

View File

@ -139,11 +139,11 @@ int main() {
|---------------------------+---------------------------------------+--------+-------------------------------------------------------------------| |---------------------------+---------------------------------------+--------+-------------------------------------------------------------------|
| ~uninitialized~ | ~int32_t~ | in | Keeps bits set for uninitialized data | | ~uninitialized~ | ~int32_t~ | in | Keeps bits set for uninitialized data |
| ~rescale_factor_ee~ | ~double~ | in | The distance scaling factor | | ~rescale_factor_ee~ | ~double~ | in | The distance scaling factor |
| ~rescale_factor_en~ | ~double[nucl_num]~ | in | The distance scaling factor | | ~rescale_factor_en~ | ~double[type_nucl_num]~ | in | The distance scaling factor |
| ~aord_num~ | ~int64_t~ | in | The number of a coeffecients | | ~aord_num~ | ~int64_t~ | in | The number of a coeffecients |
| ~bord_num~ | ~int64_t~ | in | The number of b coeffecients | | ~bord_num~ | ~int64_t~ | in | The number of b coeffecients |
| ~cord_num~ | ~int64_t~ | in | The number of c coeffecients | | ~cord_num~ | ~int64_t~ | in | The number of c coeffecients |
| ~type_nucl_num~ | ~int64_t~ | in | Number of Nucleii types | | ~type_nucl_num~ | ~int64_t~ | in | Number of Nuclei types |
| ~type_nucl_vector~ | ~int64_t[nucl_num]~ | in | IDs of types of Nuclei | | ~type_nucl_vector~ | ~int64_t[nucl_num]~ | in | IDs of types of Nuclei |
| ~a_vector~ | ~double[aord_num + 1][type_nucl_num]~ | in | a polynomial coefficients | | ~a_vector~ | ~double[aord_num + 1][type_nucl_num]~ | in | a polynomial coefficients |
| ~b_vector~ | ~double[bord_num + 1]~ | in | b polynomial coefficients | | ~b_vector~ | ~double[bord_num + 1]~ | in | b polynomial coefficients |
@ -167,7 +167,7 @@ int main() {
|-------------------------------------+-------------------------------------------------------------------+---------------------------------------------------------------------------------------------------------| |-------------------------------------+-------------------------------------------------------------------+---------------------------------------------------------------------------------------------------------|
| ~dim_c_vector~ | ~int64_t~ | Number of unique C coefficients | | ~dim_c_vector~ | ~int64_t~ | Number of unique C coefficients |
| ~dim_c_vector_date~ | ~uint64_t~ | Number of unique C coefficients | | ~dim_c_vector_date~ | ~uint64_t~ | Number of unique C coefficients |
| ~asymp_jasa~ | ~double[nucl_num]~ | Asymptotic component | | ~asymp_jasa~ | ~double[type_nucl_num]~ | Asymptotic component |
| ~asymp_jasa_date~ | ~uint64_t~ | Ladt modification of the asymptotic component | | ~asymp_jasa_date~ | ~uint64_t~ | Ladt modification of the asymptotic component |
| ~asymp_jasb~ | ~double[2]~ | Asymptotic component (up- or down-spin) | | ~asymp_jasb~ | ~double[2]~ | Asymptotic component (up- or down-spin) |
| ~asymp_jasb_date~ | ~uint64_t~ | Ladt modification of the asymptotic component | | ~asymp_jasb_date~ | ~uint64_t~ | Ladt modification of the asymptotic component |
@ -882,10 +882,10 @@ qmckl_set_jastrow_rescale_factor_en(qmckl_context context,
} }
qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero;
mem_info.size = ctx->nucleus.num * sizeof(double); mem_info.size = ctx->jastrow.type_nucl_num * sizeof(double);
ctx->jastrow.rescale_factor_en = (double*) qmckl_malloc(context, mem_info); ctx->jastrow.rescale_factor_en = (double*) qmckl_malloc(context, mem_info);
for (int64_t i=0 ; i<ctx->nucleus.num ; ++i) { for (int64_t i=0 ; i<ctx->jastrow.type_nucl_num ; ++i) {
if (rescale_factor_en[i] <= 0.0) { if (rescale_factor_en[i] <= 0.0) {
return qmckl_failwith( context, return qmckl_failwith( context,
QMCKL_INVALID_ARG_2, QMCKL_INVALID_ARG_2,
@ -948,7 +948,8 @@ qmckl_exit_code qmckl_finalize_jastrow(qmckl_context context) {
ctx->jastrow.gpu_offload = true; // ctx->electron.num > 100; ctx->jastrow.gpu_offload = true; // ctx->electron.num > 100;
#endif #endif
qmckl_exit_code rc = QMCKL_SUCCESS;
qmckl_exit_code rc = qmckl_context_touch(context);
return rc; return rc;
@ -1415,7 +1416,7 @@ qmckl_get_jastrow_rescale_factor_en (const qmckl_context context,
return QMCKL_NOT_PROVIDED; return QMCKL_NOT_PROVIDED;
} }
if (size_max < ctx->nucleus.num) { if (size_max < ctx->jastrow.type_nucl_num) {
return qmckl_failwith( context, return qmckl_failwith( context,
QMCKL_INVALID_ARG_3, QMCKL_INVALID_ARG_3,
"qmckl_get_jastrow_rescale_factor_en", "qmckl_get_jastrow_rescale_factor_en",
@ -1423,7 +1424,7 @@ qmckl_get_jastrow_rescale_factor_en (const qmckl_context context,
} }
assert(ctx->jastrow.rescale_factor_en != NULL); assert(ctx->jastrow.rescale_factor_en != NULL);
for (int64_t i=0 ; i<ctx->nucleus.num ; ++i) { for (int64_t i=0 ; i<ctx->jastrow.type_nucl_num ; ++i) {
rescale_factor_en[i] = ctx->jastrow.rescale_factor_en[i]; rescale_factor_en[i] = ctx->jastrow.rescale_factor_en[i];
} }
@ -1961,7 +1962,7 @@ assert(rc == QMCKL_SUCCESS);
double k_ee = 0.; double k_ee = 0.;
double k_en[2] = { 0., 0. }; double k_en[2] = { 0., 0. };
rc = qmckl_check(context, rc = qmckl_check(context,
qmckl_set_jastrow_rescale_factor_en(context, rescale_factor_en, nucl_num) qmckl_set_jastrow_rescale_factor_en(context, rescale_factor_en, type_nucl_num)
); );
assert(rc == QMCKL_SUCCESS); assert(rc == QMCKL_SUCCESS);
@ -1977,10 +1978,10 @@ assert(rc == QMCKL_SUCCESS);
assert(k_ee == rescale_factor_ee); assert(k_ee == rescale_factor_ee);
rc = qmckl_check(context, rc = qmckl_check(context,
qmckl_get_jastrow_rescale_factor_en (context, &(k_en[0]), nucl_num) qmckl_get_jastrow_rescale_factor_en (context, &(k_en[0]), type_nucl_num)
); );
assert(rc == QMCKL_SUCCESS); assert(rc == QMCKL_SUCCESS);
for (int i=0 ; i<nucl_num ; ++i) { for (int i=0 ; i<type_nucl_num ; ++i) {
assert(k_en[i] == rescale_factor_en[i]); assert(k_en[i] == rescale_factor_en[i]);
} }
@ -3017,7 +3018,7 @@ qmckl_get_jastrow_asymp_jasa(qmckl_context context,
qmckl_context_struct* const ctx = (qmckl_context_struct*) context; qmckl_context_struct* const ctx = (qmckl_context_struct*) context;
assert (ctx != NULL); assert (ctx != NULL);
int64_t sze = ctx->nucleus.num; int64_t sze = ctx->jastrow.type_nucl_num;
if (size_max < sze) { if (size_max < sze) {
return qmckl_failwith( context, return qmckl_failwith( context,
QMCKL_INVALID_ARG_3, QMCKL_INVALID_ARG_3,
@ -3081,7 +3082,7 @@ qmckl_exit_code qmckl_provide_jastrow_asymp_jasa(qmckl_context context)
if (ctx->jastrow.asymp_jasa == NULL) { if (ctx->jastrow.asymp_jasa == NULL) {
qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero;
mem_info.size = ctx->nucleus.num * sizeof(double); mem_info.size = ctx->jastrow.type_nucl_num * sizeof(double);
double* asymp_jasa = (double*) qmckl_malloc(context, mem_info); double* asymp_jasa = (double*) qmckl_malloc(context, mem_info);
if (asymp_jasa == NULL) { if (asymp_jasa == NULL) {
@ -3097,7 +3098,7 @@ qmckl_exit_code qmckl_provide_jastrow_asymp_jasa(qmckl_context context)
ctx->jastrow.aord_num, ctx->jastrow.aord_num,
ctx->jastrow.a_vector, ctx->jastrow.a_vector,
ctx->jastrow.rescale_factor_en, ctx->jastrow.rescale_factor_en,
ctx->nucleus.num, ctx->jastrow.type_nucl_num,
ctx->jastrow.asymp_jasa); ctx->jastrow.asymp_jasa);
if (rc != QMCKL_SUCCESS) { if (rc != QMCKL_SUCCESS) {
return rc; return rc;
@ -3119,26 +3120,26 @@ qmckl_exit_code qmckl_provide_jastrow_asymp_jasa(qmckl_context context)
#+NAME: qmckl_asymp_jasa_args #+NAME: qmckl_asymp_jasa_args
| Variable | Type | In/Out | Description | | Variable | Type | In/Out | Description |
|---------------------+----------------------+--------+----------------------------| |---------------------+-------------------------------------+--------+----------------------------|
| ~context~ | ~qmckl_context~ | in | Global state | | ~context~ | ~qmckl_context~ | in | Global state |
| ~aord_num~ | ~int64_t~ | in | Order of the polynomial | | ~aord_num~ | ~int64_t~ | in | Order of the polynomial |
| ~a_vector~ | ~double[aord_num+1]~ | in | Values of a | | ~a_vector~ | ~double[type_nucl_num][aord_num+1]~ | in | Values of a |
| ~rescale_factor_en~ | ~double~ | in | Electron nucleus distances | | ~rescale_factor_en~ | ~double~ | in | Electron nucleus distances |
| ~nucl_num~ | ~int64_t~ | in | Number of nuclei | | ~type_nucl_num~ | ~int64_t~ | in | Number of nucleus types |
| ~asymp_jasa~ | ~double[nucl_num]~ | out | Asymptotic value | | ~asymp_jasa~ | ~double[type_nucl_num]~ | out | Asymptotic value |
#+begin_src f90 :comments org :tangle (eval f) :noweb yes #+begin_src f90 :comments org :tangle (eval f) :noweb yes
integer function qmckl_compute_jastrow_asymp_jasa_f(context, aord_num, a_vector, & integer function qmckl_compute_jastrow_asymp_jasa_f(context, aord_num, a_vector, &
rescale_factor_en, nucl_num, asymp_jasa) & rescale_factor_en, type_nucl_num, asymp_jasa) &
result(info) result(info)
use qmckl use qmckl
implicit none implicit none
integer(qmckl_context), intent(in) :: context integer(qmckl_context), intent(in) :: context
integer*8 , intent(in) :: aord_num integer*8 , intent(in) :: aord_num
double precision , intent(in) :: a_vector(aord_num + 1) integer*8 , intent(in) :: type_nucl_num
integer*8 , intent(in) :: nucl_num double precision , intent(in) :: a_vector(aord_num + 1, type_nucl_num)
double precision , intent(in) :: rescale_factor_en(nucl_num) double precision , intent(in) :: rescale_factor_en(type_nucl_num)
double precision , intent(out) :: asymp_jasa(nucl_num) double precision , intent(out) :: asymp_jasa(type_nucl_num)
integer*8 :: i, j, p integer*8 :: i, j, p
double precision :: kappa_inv, x, asym_one double precision :: kappa_inv, x, asym_one
@ -3156,16 +3157,16 @@ integer function qmckl_compute_jastrow_asymp_jasa_f(context, aord_num, a_vector,
return return
endif endif
do i=1,nucl_num do i=1,type_nucl_num
kappa_inv = 1.0d0 / rescale_factor_en(i) kappa_inv = 1.0d0 / rescale_factor_en(i)
asymp_jasa(i) = a_vector(1) * kappa_inv / (1.0d0 + a_vector(2) * kappa_inv) asymp_jasa(i) = a_vector(1,i) * kappa_inv / (1.0d0 + a_vector(2,i) * kappa_inv)
x = kappa_inv x = kappa_inv
do p = 1, aord_num do p = 1, aord_num
x = x * kappa_inv x = x * kappa_inv
asymp_jasa(i) = asymp_jasa(i) + a_vector(p + 1) * x asymp_jasa(i) = asymp_jasa(i) + a_vector(p + 1, i) * x
end do end do
end do end do
@ -3179,7 +3180,7 @@ qmckl_exit_code qmckl_compute_jastrow_asymp_jasa (
const int64_t aord_num, const int64_t aord_num,
const double* a_vector, const double* a_vector,
double* const rescale_factor_en, double* const rescale_factor_en,
const int64_t nucl_num, const int64_t type_nucl_num,
double* const asymp_jasa ) { double* const asymp_jasa ) {
if (context == QMCKL_NULL_CONTEXT){ if (context == QMCKL_NULL_CONTEXT){
@ -3190,14 +3191,14 @@ qmckl_exit_code qmckl_compute_jastrow_asymp_jasa (
return QMCKL_INVALID_ARG_2; return QMCKL_INVALID_ARG_2;
} }
for (int i = 0 ; i <= nucl_num; ++i) { for (int i = 0 ; i <= type_nucl_num; ++i) {
const double kappa_inv = 1.0 / rescale_factor_en[i]; const double kappa_inv = 1.0 / rescale_factor_en[i];
asymp_jasa[i] = a_vector[0] * kappa_inv / (1.0 + a_vector[1] * kappa_inv); asymp_jasa[i] = a_vector[0 + aord_num*i] * kappa_inv / (1.0 + a_vector[1 + aord_num*i] * kappa_inv);
double x = kappa_inv; double x = kappa_inv;
for (int p = 1; p < aord_num; ++p){ for (int p = 1; p < aord_num; ++p){
x *= kappa_inv; x *= kappa_inv;
asymp_jasa[i] = asymp_jasa[i] + a_vector[p + 1] * x; asymp_jasa[i] = asymp_jasa[i] + a_vector[p + 1 + aord_num*i] * x;
} }
} }
@ -3213,7 +3214,7 @@ qmckl_exit_code qmckl_compute_jastrow_asymp_jasa (
const int64_t aord_num, const int64_t aord_num,
const double* a_vector, const double* a_vector,
double* const rescale_factor_en, double* const rescale_factor_en,
const int64_t nucl_num, const int64_t type_nucl_num,
double* const asymp_jasa ); double* const asymp_jasa );
#+end_src #+end_src
@ -3556,14 +3557,13 @@ qmckl_exit_code qmckl_compute_factor_en (
factor_en[nw] = 0.0; factor_en[nw] = 0.0;
for (int a = 0; a < nucl_num; ++a ) { for (int a = 0; a < nucl_num; ++a ) {
for (int i = 0; i < elec_num; ++i ) { for (int i = 0; i < elec_num; ++i ) {
// x = ee_distance_rescaled[j * (walk_num * elec_num) + i * (walk_num) + nw];
x = en_distance_rescaled[i + a * elec_num + nw * (elec_num * nucl_num)]; x = en_distance_rescaled[i + a * elec_num + nw * (elec_num * nucl_num)];
x1 = x; x1 = x;
power_ser = 0.0; power_ser = 0.0;
for (int p = 2; p < aord_num+1; ++p) { for (int p = 2; p < aord_num+1; ++p) {
x = x * x1; x = x * x1;
power_ser = power_ser + a_vector[(p+1)-1 + (type_nucl_vector[a]-1) * aord_num] * x; power_ser = power_ser + a_vector[p+ (type_nucl_vector[a]-1) * aord_num] * x;
} }
factor_en[nw] = factor_en[nw] + a_vector[0 + (type_nucl_vector[a]-1)*aord_num] * x1 / \ factor_en[nw] = factor_en[nw] + a_vector[0 + (type_nucl_vector[a]-1)*aord_num] * x1 / \
@ -5576,6 +5576,8 @@ qmckl_exit_code qmckl_provide_en_distance_rescaled(qmckl_context context)
qmckl_compute_en_distance_rescaled(context, qmckl_compute_en_distance_rescaled(context,
ctx->electron.num, ctx->electron.num,
ctx->nucleus.num, ctx->nucleus.num,
ctx->jastrow.type_nucl_num,
ctx->jastrow.type_nucl_vector,
ctx->jastrow.rescale_factor_en, ctx->jastrow.rescale_factor_en,
ctx->electron.walker.num, ctx->electron.walker.num,
ctx->electron.walker.point.coord.data, ctx->electron.walker.point.coord.data,
@ -5605,15 +5607,17 @@ qmckl_exit_code qmckl_provide_en_distance_rescaled(qmckl_context context)
| ~context~ | ~qmckl_context~ | in | Global state | | ~context~ | ~qmckl_context~ | in | Global state |
| ~elec_num~ | ~int64_t~ | in | Number of electrons | | ~elec_num~ | ~int64_t~ | in | Number of electrons |
| ~nucl_num~ | ~int64_t~ | in | Number of nuclei | | ~nucl_num~ | ~int64_t~ | in | Number of nuclei |
| ~rescale_factor_en~ | ~double[nucl_num]~ | in | The factor for rescaled distances | | ~type_nucl_num~ | ~int64_t~ | in | Number of types of nuclei |
| ~type_nucl_vector~ | ~int64_t[nucl_num]~ | in | Number of types of nuclei |
| ~rescale_factor_en~ | ~double[type_nucl_num]~ | in | The factor for rescaled distances |
| ~walk_num~ | ~int64_t~ | in | Number of walkers | | ~walk_num~ | ~int64_t~ | in | Number of walkers |
| ~elec_coord~ | ~double[3][walk_num][elec_num]~ | in | Electron coordinates | | ~elec_coord~ | ~double[3][walk_num][elec_num]~ | in | Electron coordinates |
| ~nucl_coord~ | ~double[3][elec_num]~ | in | Nuclear coordinates | | ~nucl_coord~ | ~double[3][elec_num]~ | in | Nuclear coordinates |
| ~en_distance_rescaled~ | ~double[walk_num][nucl_num][elec_num]~ | out | Electron-nucleus distances | | ~en_distance_rescaled~ | ~double[walk_num][nucl_num][elec_num]~ | out | Electron-nucleus distances |
| | | | |
#+begin_src f90 :comments org :tangle (eval f) :noweb yes #+begin_src f90 :comments org :tangle (eval f) :noweb yes
integer function qmckl_compute_en_distance_rescaled_f(context, elec_num, nucl_num, rescale_factor_en, walk_num, elec_coord, & integer function qmckl_compute_en_distance_rescaled_f(context, elec_num, nucl_num, type_nucl_num, &
type_nucl_vector, rescale_factor_en, walk_num, elec_coord, &
nucl_coord, en_distance_rescaled) & nucl_coord, en_distance_rescaled) &
result(info) result(info)
use qmckl use qmckl
@ -5621,7 +5625,9 @@ integer function qmckl_compute_en_distance_rescaled_f(context, elec_num, nucl_nu
integer(qmckl_context), intent(in) :: context integer(qmckl_context), intent(in) :: context
integer*8 , intent(in) :: elec_num integer*8 , intent(in) :: elec_num
integer*8 , intent(in) :: nucl_num integer*8 , intent(in) :: nucl_num
double precision , intent(in) :: rescale_factor_en(nucl_num) integer*8 , intent(in) :: type_nucl_num
integer*8 , intent(in) :: type_nucl_vector(nucl_num)
double precision , intent(in) :: rescale_factor_en(type_nucl_num)
integer*8 , intent(in) :: walk_num integer*8 , intent(in) :: walk_num
double precision , intent(in) :: elec_coord(elec_num,walk_num,3) double precision , intent(in) :: elec_coord(elec_num,walk_num,3)
double precision , intent(in) :: nucl_coord(nucl_num,3) double precision , intent(in) :: nucl_coord(nucl_num,3)
@ -5657,7 +5663,7 @@ integer function qmckl_compute_en_distance_rescaled_f(context, elec_num, nucl_nu
do k=1,walk_num do k=1,walk_num
info = qmckl_distance_rescaled(context, 'T', 'T', elec_num, 1_8, & info = qmckl_distance_rescaled(context, 'T', 'T', elec_num, 1_8, &
elec_coord(1,k,1), elec_num*walk_num, coord, 1_8, & elec_coord(1,k,1), elec_num*walk_num, coord, 1_8, &
en_distance_rescaled(1,i,k), elec_num, rescale_factor_en(i)) en_distance_rescaled(1,i,k), elec_num, rescale_factor_en(type_nucl_vector(i)))
if (info /= QMCKL_SUCCESS) then if (info /= QMCKL_SUCCESS) then
return return
endif endif
@ -5672,6 +5678,8 @@ qmckl_exit_code qmckl_compute_en_distance_rescaled (
const qmckl_context context, const qmckl_context context,
const int64_t elec_num, const int64_t elec_num,
const int64_t nucl_num, const int64_t nucl_num,
const int64_t type_nucl_num,
int64_t* const type_nucl_vector,
const double* rescale_factor_en, const double* rescale_factor_en,
const int64_t walk_num, const int64_t walk_num,
const double* elec_coord, const double* elec_coord,
@ -5687,6 +5695,8 @@ qmckl_exit_code qmckl_compute_en_distance_rescaled (
(context, & (context, &
elec_num, & elec_num, &
nucl_num, & nucl_num, &
type_nucl_num, &
type_nucl_vector, &
rescale_factor_en, & rescale_factor_en, &
walk_num, & walk_num, &
elec_coord, & elec_coord, &
@ -5700,7 +5710,9 @@ qmckl_exit_code qmckl_compute_en_distance_rescaled (
integer (c_int64_t) , intent(in) , value :: context integer (c_int64_t) , intent(in) , value :: context
integer (c_int64_t) , intent(in) , value :: elec_num integer (c_int64_t) , intent(in) , value :: elec_num
integer (c_int64_t) , intent(in) , value :: nucl_num integer (c_int64_t) , intent(in) , value :: nucl_num
real (c_double ) , intent(in) :: rescale_factor_en(nucl_num) integer (c_int64_t) , intent(in) , value :: type_nucl_num
integer (c_int64_t) , intent(in) :: type_nucl_vector(nucl_num)
real (c_double ) , intent(in) :: rescale_factor_en(type_nucl_num)
integer (c_int64_t) , intent(in) , value :: walk_num integer (c_int64_t) , intent(in) , value :: walk_num
real (c_double ) , intent(in) :: elec_coord(elec_num,walk_num,3) real (c_double ) , intent(in) :: elec_coord(elec_num,walk_num,3)
real (c_double ) , intent(in) :: nucl_coord(elec_num,3) real (c_double ) , intent(in) :: nucl_coord(elec_num,3)
@ -5711,6 +5723,8 @@ qmckl_exit_code qmckl_compute_en_distance_rescaled (
(context, & (context, &
elec_num, & elec_num, &
nucl_num, & nucl_num, &
type_nucl_num, &
type_nucl_vector, &
rescale_factor_en, & rescale_factor_en, &
walk_num, & walk_num, &
elec_coord, & elec_coord, &

View File

@ -297,8 +297,7 @@ qmckl_exit_code qmckl_finalize_mo_basis(qmckl_context context) {
ctx->mo_basis.mo_value = NULL; ctx->mo_basis.mo_value = NULL;
ctx->mo_basis.mo_value_date = 0; ctx->mo_basis.mo_value_date = 0;
} }
return qmckl_context_touch(context);
return QMCKL_SUCCESS;
} }
#+end_src #+end_src