Form of the functional
W-SLDA codes minimize functional of generic form:
where:
- is energy density functional which defines the physical system,
- is spin dependent external potential, and are chemical potentials (Lagrange multipliers) for constraining particle number.
- is external pairing potential,
- is external velocity field.
W-SLDA toolkit provides flexible framework that allows for custom definition of all these terms. User must provide body of following functions contained in file: problem-definition.h. Each function can be parametrized by User defined parameters.
Definition of the external potential
/**
* EXTERNAL POTENTIAL V_ext
* @param ix x-coordinate from range [0,NX), to convert to Cartesian use: x = DX*(ix-NX/2)
* @param iy y-coordinate from range [0,NY), to convert to Cartesian use: y = DY*(iy-NY/2),
* NOTE: in case of 1d code iy=0
* @param iz z-coordinate from range [0,NZ), to convert to Cartesian use: z = DZ*(iz-NZ/2)
* NOTE: in case of 1d and 2d codes iz=0
* @param it iteration number
* @param spin spin indicator, value from set {SPINA,SPINB}
* @param params array of input parameters, before call of this routine the params array is processed by process_params() routine
* @param extra_data_size size of extra_data in bytes, if extra_data size=0 the optional data is not uploaded
* @param extra_data optional set of data uploaded by load_extra_data()
* @return value of the external potential V_spin(x,y,z)
* */
double v_ext(int ix, int iy, int iz, int it, int spin, double *params, size_t extra_data_size, void *extra_data)
{
// ADD HERE FORMULA FOR V_ext(r)
double V_ext = 0.0;
return V_ext;
}
Definition of the external pairing potential
/**
* EXTERNAL PAIRING POTENTIAL Delta_ext
* @param ix x-coordinate from range [0,NX), to convert to Cartesian use: x = DX*(ix-NX/2)
* @param iy y-coordinate from range [0,NY), to convert to Cartesian use: y = DY*(iy-NY/2)
* NOTE: in case of 1d code iy=0
* @param iz z-coordinate from range [0,NZ), to convert to Cartesian use: z = DZ*(iz-NZ/2)
* NOTE: in case of 1d and 2d codes iz=0
* @param it iteration number
* @param delta - value of delta computed self-consistently for given iteration it.
* @param params array of input parameters, before call of this routine the params array is processed by process_params() routine
* @param extra_data_size size of extra_data in bytes, if extra_data size=0 the optional data is not uploaded
* @param extra_data optional set of data uploaded by load_extra_data()
* @return value of external pairing potential Delta_{ext}(x,y,z)
* */
double complex delta_ext(int ix, int iy, int iz, int it, double complex delta, double *params, size_t extra_data_size, void *extra_data)
{
// ADD HERE FORMULA FOR Delta_ext(r)
double complex D_ext = 0.0 + I*0.0;
return D_ext;
}
Note: Due to technical reasons this function differs with respect to return type between st-wslda
and td-wslda
codes. Namly:
-
st-wslda
: return type must be C99 double complex -
st-wslda
: return type must be compatible with CUDA Complex
Definition of the external velocity field
/**
* EXTERNAL VELOCITY FIELD vec[v]_ext
* @param ix x-coordinate from range [0,NX), to convert to Cartesian use: x = DX*(ix-NX/2)
* @param iy y-coordinate from range [0,NY), to convert to Cartesian use: y = DY*(iy-NY/2)
* NOTE: in case of 1d code iy=0
* @param iz z-coordinate from range [0,NZ), to convert to Cartesian use: z = DZ*(iz-NZ/2)
* NOTE: in case of 1d and 2d codes iz=0
* @param it iteration number
* @param spin spin indicator, value from set {SPINA,SPINB}
* @param coordinate - Cartesian coordinate of the external velocity vector that should be computed, value from set {XAXIS, YAXIS, ZAXIS}
* NOTE: for 1d code only XAXIS is requested, for 2d code XAXIS and YAXIS are requested.
* @param params array of input parameters, before call of this routine the params array is processed by process_params() routine
* @param extra_data_size size of extra_data in bytes, if extra_data size=0 the optional data is not uploaded
* @param extra_data optional set of data uploaded by load_extra_data()
* @return value of the external velocity vector v_ext(x,y,z)
* */
double velocity_ext(int ix, int iy, int iz, int it, int spin, int coordinate, double *params, size_t extra_data_size, void *extra_data)
{
// ADD HERE FORMULAS FOR vec{v}_ext=(vx, vy, vz)
double v_ext;
if(coordinate==XAXIS) v_ext=0.0;
if(coordinate==YAXIS) v_ext=0.0;
if(coordinate==ZAXIS) v_ext=0.0;
return v_ext;
}