Represents an MBD calculation.
Initialize an MBD calculation from an MBD input.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(inout), | target | :: | this | ||
type(mbd_input_t), | intent(in) | :: | input | MBD input. |
Finalize an MBD calculation.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(inout), | target | :: | this |
Update whether to calculate forces.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(inout) | :: | this | |||
logical, | intent(in) | :: | forces | Whether to calcualte forces. |
Update atomic coordinates.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(inout) | :: | this | |||
real(kind=dp), | intent(in) | :: | coords(:,:) | (, a.u.) New atomic coordinates. |
Update unit-cell lattice vectors.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(inout) | :: | this | |||
real(kind=dp), | intent(in) | :: | latt_vecs(:,:) | (, a.u.) New lattice vectors in columns. |
Update vdW parameters in a custom way.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(inout) | :: | this | |||
real(kind=dp), | intent(in) | :: | alpha_0(:) | (a.u.) New atomic static polarizabilities. |
||
real(kind=dp), | intent(in) | :: | C6(:) | (a.u.) New atomic coefficients. |
||
real(kind=dp), | intent(in) | :: | r_vdw(:) | (a.u.) New atomic vdW radii. |
Update vdW parameters based on scaling of free-atom values.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(inout) | :: | this | |||
real(kind=dp), | intent(in) | :: | ratios(:) | Ratios of atomic volumes in the system and in vacuum. |
Update vdW parameters for the MBD-NL method.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(inout) | :: | this | |||
real(kind=dp), | intent(in) | :: | alpha_0_ratios(:) | Ratios of free-atom exact static polarizabilities and those from the VV functional. |
||
real(kind=dp), | intent(in) | :: | C6_ratios(:) | Ratios of free-atom exact coefficients and those from the VV functional. |
Evaluate a given vdW method for a given system and vdW parameters, retrieve energy.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(inout) | :: | this | |||
real(kind=dp), | intent(out) | :: | energy | (a.u.) VdW energy. |
Retrieve nuclear energy gradients if they were requested in the MBD input.
The gradients are calculated together with the energy, so a call to this method must be preceeded by a call to evaluate_vdw_method. For the same reason, the gradients must be requested prior to this called via calculate_forces.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(in) | :: | this | |||
real(kind=dp), | intent(out) | :: | gradients(:,:) | (, a.u.) Energy gradients, , index runs over columns. |
Get gradients of the energy w.r.t. Hirshfeld ratios if they were requested in the MBD input.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(inout) | :: | this | |||
real(kind=dp), | intent(out) | :: | dE_dratios(:) | Gradients of the energy w.r.t. Hirshfeld ratios. |
Provide lattice-vector energy gradients if they were requested in the MBD input.
The gradients are actually calculated together with the energy, so a call to this method must be preceeded by a call to evaluate_vdw_method. For the same reason, the gradients must be requested prior to this called via calculate_forces.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(in) | :: | this | |||
real(kind=dp), | intent(out) | :: | latt_derivs(:,:) | (, a.u.) Energy gradients, , index runs over columns. |
Provide stress tensor of the lattice.
This is a utility function wrapping get_lattice_derivs. The lattice vector gradients are coverted to the stress tensor.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(in) | :: | this | |||
real(kind=dp), | intent(out) | :: | stress(:,:) | (, a.u.) Stress tensor. |
Provide MBD spectrum if it was requested in the MBD input.
The spectrum is actually calculated together with the energy, so a call to this method must be preceeded by a call to evaluate_vdw_method. For the same reason, the spectrum must be requested prior to this call via calculate_spectrum.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(inout) | :: | this | |||
real(kind=dp), | intent(out) | :: | spectrum(:) | (, a.u.) Energies (frequencies) of coupled MBD modues, . |
||
real(kind=dp), | intent(out), | optional | allocatable | :: | modes(:,:) | () Coupled-mode wave functions (MBD eigenstates), , in the basis of uncoupled states, , index runs over columns. To save memory, the argument must be allocatable, and the method transfers allocation from the internal state to the argument. For this reason, the method can be called only once wih this optional argument per calculation. |
Provide RPA orders if they were requested in the MBD input.
The orders are actually calculated together with the energy, so a call to this method must be preceeded by a call to evaluate_vdw_method. For the same reason, the spectrum must be requested prior to this call via do_rpa and rpa_orders.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(inout) | :: | this | |||
real(kind=dp), | intent(out), | allocatable | :: | rpa_orders(:) | (a.u.) MBD energy decomposed to RPA orders. |
Retrieve an exception in the MBD calculation if it occured.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(mbd_calc_t), | intent(inout) | :: | this | |||
integer, | intent(out) | :: | code | Exception code, values defined in mbd_constants. |
||
character(len=*), | intent(out) | :: | origin | Exception origin. |
||
character(len=*), | intent(out) | :: | msg | Exception message. |