Octopus
mix_oct_m Module Reference

Data Types

type  mix_t
 God class for mixing. More...
 
interface  mixfield_get_vnew
 
type  mixfield_ptr_t
 
interface  mixfield_set_vin
 
interface  mixfield_set_vout
 
type  mixfield_t
 Quantities used in mixing: Input, output and new potentials, and the residuals. More...
 

Functions/Subroutines

subroutine, public mix_init (smix, namespace, space, der, d1, d2, def_, func_type_, prefix_)
 Initialise mix_t instance. More...
 
subroutine, public mix_clear (smix)
 
subroutine, public mix_end (smix)
 
subroutine, public mix_dump (namespace, restart, smix, space, mesh, ierr)
 
subroutine, public mix_load (namespace, restart, smix, space, mesh, ierr)
 
real(real64) pure function, public mix_coefficient (this)
 
integer pure function, public mix_scheme (this)
 
integer pure function, public mix_d3 (this)
 
subroutine, public mix_get_field (this, mixfield)
 
subroutine, public mixing (namespace, smix)
 Main entry-point to SCF mixer. More...
 
subroutine, public mix_add_auxmixfield (namespace, smix, mixfield)
 
subroutine, public mixfield_init (smix, mixfield, d1, d2, d3, func_type)
 Initialise all attributes of a mixfield instance. More...
 
subroutine, public mixfield_end (smix, mixfield)
 Deallocate all arrays of a mixfield instance. More...
 
subroutine, public mixfield_clear (scheme, mixfield)
 Zero all potential and field attributes of a mixfield instance. More...
 
subroutine ddmixfield_set_vin (mixfield, vin_re, vin_im)
 
subroutine ddmixfield_set_vout (mixfield, vout_re, vout_im)
 
subroutine ddmixfield_get_vnew (mixfield, re, im)
 
subroutine compute_residuals_aux_field (this)
 Compute the residuals of the auxilliary fields. More...
 
subroutine linear_mixing_aux_field (smix)
 Linear mixing of the auxilliary fields. More...
 
subroutine, public dmixing (namespace, smix, vin, vout, vnew)
 Mix the input and output potentials (or densities) according to some scheme. More...
 
subroutine dmixing_linear (d1, d2, coeff, vin, residual, vnew)
 Linear mixing of the input and output potentials. More...
 
subroutine dcompute_residual (smix, dims, vin, vout, residual)
 Compute the residual of the potential (or density) for SCF mixing. More...
 
subroutine dmixing_broyden (namespace, smix, f, vin, vnew, iter)
 Top-level routine for Broyden mixing. More...
 
subroutine dbroyden_extrapolation (this, namespace, coeff, d1, d2, vin, vnew, iter_used, f, df, dv)
 Broyden mixing implementation. More...
 
subroutine dbroyden_extrapolation_aux (this, namespace, ii, coeff, iter_used, dbeta, dwork, zbeta, zwork)
 Broyden mixing for auxilliary fields. More...
 
subroutine dmixing_diis (this, vin, residual, vnew, iter)
 Direct inversion in the iterative subspace. More...
 
subroutine dmixing_build_matrix (this, df, size, d2, ww, beta)
 
real(real64) function dmix_dotp (this, xx, yy, reduce)
 
subroutine dmixfield_set_vin (mixfield, vin)
 
subroutine dmixfield_set_vout (mixfield, vout)
 
subroutine dmixfield_get_vnew (mixfield, vnew)
 
subroutine, public zmixing (namespace, smix, vin, vout, vnew)
 Mix the input and output potentials (or densities) according to some scheme. More...
 
subroutine zmixing_linear (d1, d2, coeff, vin, residual, vnew)
 Linear mixing of the input and output potentials. More...
 
subroutine zcompute_residual (smix, dims, vin, vout, residual)
 Compute the residual of the potential (or density) for SCF mixing. More...
 
subroutine zmixing_broyden (namespace, smix, f, vin, vnew, iter)
 Top-level routine for Broyden mixing. More...
 
subroutine zbroyden_extrapolation (this, namespace, coeff, d1, d2, vin, vnew, iter_used, f, df, dv)
 Broyden mixing implementation. More...
 
subroutine zbroyden_extrapolation_aux (this, namespace, ii, coeff, iter_used, dbeta, dwork, zbeta, zwork)
 Broyden mixing for auxilliary fields. More...
 
subroutine zmixing_diis (this, vin, residual, vnew, iter)
 Direct inversion in the iterative subspace. More...
 
subroutine zmixing_build_matrix (this, df, size, d2, ww, beta)
 
complex(real64) function zmix_dotp (this, xx, yy, reduce)
 
subroutine zmixfield_set_vin (mixfield, vin)
 
subroutine zmixfield_set_vout (mixfield, vout)
 
subroutine zmixfield_get_vnew (mixfield, vnew)
 

Variables

integer, parameter max_auxmixfield = 5
 

Function/Subroutine Documentation

◆ mix_init()

subroutine, public mix_oct_m::mix_init ( type(mix_t), intent(out)  smix,
type(namespace_t), intent(in)  namespace,
class(space_t), intent(in)  space,
type(derivatives_t), intent(in), target  der,
integer, intent(in)  d1,
integer, intent(in)  d2,
integer, intent(in), optional  def_,
type(type_t), intent(in), optional  func_type_,
character(len=*), intent(in), optional  prefix_ 
)

Initialise mix_t instance.

Parameters
[in]d1Limit of the potential, typically meshnp
[in]d2Limits of the potential, typically stdnspin

Definition at line 264 of file mix.F90.

◆ mix_clear()

subroutine, public mix_oct_m::mix_clear ( type(mix_t), intent(inout)  smix)

Definition at line 506 of file mix.F90.

◆ mix_end()

subroutine, public mix_oct_m::mix_end ( type(mix_t), intent(inout)  smix)

Definition at line 521 of file mix.F90.

◆ mix_dump()

subroutine, public mix_oct_m::mix_dump ( type(namespace_t), intent(in)  namespace,
type(restart_t), intent(in)  restart,
type(mix_t), intent(in)  smix,
class(space_t), intent(in)  space,
class(mesh_t), intent(in)  mesh,
integer, intent(out)  ierr 
)

Definition at line 542 of file mix.F90.

◆ mix_load()

subroutine, public mix_oct_m::mix_load ( type(namespace_t), intent(in)  namespace,
type(restart_t), intent(in)  restart,
type(mix_t), intent(inout)  smix,
class(space_t), intent(in)  space,
class(mesh_t), intent(in)  mesh,
integer, intent(out)  ierr 
)

Definition at line 644 of file mix.F90.

◆ mix_coefficient()

real(real64) pure function, public mix_oct_m::mix_coefficient ( type(mix_t), intent(in)  this)

Definition at line 773 of file mix.F90.

◆ mix_scheme()

integer pure function, public mix_oct_m::mix_scheme ( type(mix_t), intent(in)  this)

Definition at line 779 of file mix.F90.

◆ mix_d3()

integer pure function, public mix_oct_m::mix_d3 ( type(mix_t), intent(in)  this)

Definition at line 785 of file mix.F90.

◆ mix_get_field()

subroutine, public mix_oct_m::mix_get_field ( type(mix_t), intent(in), target  this,
type(mixfield_t), intent(out), pointer  mixfield 
)

Definition at line 791 of file mix.F90.

◆ mixing()

subroutine, public mix_oct_m::mixing ( type(namespace_t), intent(in)  namespace,
type(mix_t), intent(inout)  smix 
)

Main entry-point to SCF mixer.

Definition at line 799 of file mix.F90.

◆ mix_add_auxmixfield()

subroutine, public mix_oct_m::mix_add_auxmixfield ( type(namespace_t), intent(in)  namespace,
type(mix_t), intent(inout)  smix,
type(mixfield_t), intent(in), target  mixfield 
)

Definition at line 814 of file mix.F90.

◆ mixfield_init()

subroutine, public mix_oct_m::mixfield_init ( type(mix_t), intent(in)  smix,
type(mixfield_t), intent(inout)  mixfield,
integer, intent(in)  d1,
integer, intent(in)  d2,
integer, intent(in)  d3,
type(type_t), intent(in)  func_type 
)

Initialise all attributes of a mixfield instance.

Definition at line 833 of file mix.F90.

◆ mixfield_end()

subroutine, public mix_oct_m::mixfield_end ( type(mix_t), intent(inout)  smix,
type(mixfield_t), intent(inout)  mixfield 
)

Deallocate all arrays of a mixfield instance.

Definition at line 877 of file mix.F90.

◆ mixfield_clear()

subroutine, public mix_oct_m::mixfield_clear ( integer, intent(in)  scheme,
type(mixfield_t), intent(inout)  mixfield 
)

Zero all potential and field attributes of a mixfield instance.

Definition at line 914 of file mix.F90.

◆ ddmixfield_set_vin()

subroutine mix_oct_m::ddmixfield_set_vin ( type(mixfield_t), intent(inout)  mixfield,
real(real64), dimension(:,:), intent(in)  vin_re,
real(real64), dimension(:,:), intent(in)  vin_im 
)
private

Definition at line 960 of file mix.F90.

◆ ddmixfield_set_vout()

subroutine mix_oct_m::ddmixfield_set_vout ( type(mixfield_t), intent(inout)  mixfield,
real(real64), dimension(:,:), intent(in)  vout_re,
real(real64), dimension(:,:), intent(in)  vout_im 
)
private

Definition at line 972 of file mix.F90.

◆ ddmixfield_get_vnew()

subroutine mix_oct_m::ddmixfield_get_vnew ( type(mixfield_t), intent(in)  mixfield,
real(real64), dimension(:,:), intent(inout)  re,
real(real64), dimension(:,:), intent(inout)  im 
)
private

Definition at line 984 of file mix.F90.

◆ compute_residuals_aux_field()

subroutine mix_oct_m::compute_residuals_aux_field ( class(mix_t), intent(inout)  this)
private

Compute the residuals of the auxilliary fields.

Given potential (or density) attributes of an smixauxmixfield(i) instance, compute the residuals.

The if-statement is required because array attributes of an auxmixfield instance can be contain either real or complex. Whether the real or complex arrays should be used is indicated by the attribute smixauxmixfield(i)pfunc_type.

Parameters
[in,out]thisMixing instance and derivatives

Definition at line 1004 of file mix.F90.

◆ linear_mixing_aux_field()

subroutine mix_oct_m::linear_mixing_aux_field ( type(mix_t), intent(inout)  smix)
private

Linear mixing of the auxilliary fields.

Given potential (or density) attributes of an smixauxmixfield(i) instance, compute the residuals, then linearly mix the potential (or density).

The if-statement is required because array attributes of an auxmixfield instance can be contain either real or complex. Whether the real or complex arrays should be used is indicated by the attribute smixauxmixfield(i)pfunc_type.

Note
This could be changed to a type-bound method as it mutates smixauxmixfield(i)p attributes
Parameters
[in,out]smixMixing instance and derivatives

Definition at line 1041 of file mix.F90.

◆ dmixing()

subroutine, public mix_oct_m::dmixing ( type(namespace_t), intent(in)  namespace,
type(mix_t), intent(inout), target  smix,
real(real64), dimension(:, :), intent(in), contiguous  vin,
real(real64), dimension(:, :), intent(in), contiguous  vout,
real(real64), dimension(:, :), intent(out), contiguous  vnew 
)

Mix the input and output potentials (or densities) according to some scheme.

Whether \(v\) is the potential or the density is determined at run time.

Parameters
[in,out]smixMixing settings and derivatives
[in]vinInput potential
[in]voutOutput potential
[out]vnewUpdated potential

Definition at line 1140 of file mix.F90.

◆ dmixing_linear()

subroutine mix_oct_m::dmixing_linear ( integer, intent(in)  d1,
integer, intent(in)  d2,
real(real64), intent(in)  coeff,
real(real64), dimension(:, :), intent(in)  vin,
real(real64), dimension(:, :), intent(in)  residual,
real(real64), dimension(:, :), intent(out)  vnew 
)
private

Linear mixing of the input and output potentials.

Linear mixing is defined as:

\[ V_{\text{new}} = (1 - a) V_{\text{in}} + a V_{\text{out}}. \]

Using the definition of the residual, f$ V_{\text{out}} = V_{\text{in}} + R

Parameters
[in]d2Limits for each dimension of the potential
[in]coeffMixing coefficient
[in]vinInput potential
[out]vnewUpdated potential

Definition at line 1194 of file mix.F90.

◆ dcompute_residual()

subroutine mix_oct_m::dcompute_residual ( type(mix_t), intent(in)  smix,
integer, dimension(2), intent(in)  dims,
real(real64), dimension(:, :), intent(in)  vin,
real(real64), dimension(:, :), intent(in)  vout,
real(real64), dimension(:, :), intent(out), contiguous  residual 
)
private

Compute the residual of the potential (or density) for SCF mixing.

In the absence of any preconditioning, the residual is defined as \( V_\text{out} - V_\text{in}\). If preconditioning is used, please refer to the specific scheme for the definition.

Note, d1, d2 are passed in because this routine should be compatible with both the potential and and the auxiliary fields (which have differing sizes).

Parameters
[in]smixMixing settings and derivatives
[in]dimsDimensions of the terms to mix
[in]vinInput potential
[in]voutOutput potential

Definition at line 1226 of file mix.F90.

◆ dmixing_broyden()

subroutine mix_oct_m::dmixing_broyden ( type(namespace_t), intent(in)  namespace,
type(mix_t), intent(inout)  smix,
real(real64), dimension(:, :), intent(in), contiguous  f,
real(real64), dimension(:, :), intent(in), contiguous  vin,
real(real64), dimension(:, :), intent(out), contiguous  vnew,
integer, intent(in)  iter 
)
private

Top-level routine for Broyden mixing.

Performs updating of smixmixfield attributes. Whether \(v\) is the potential or the density is determined at run time.

Parameters
[in,out]smixMixing settings and derivatives
[in]fResidual
[in]vinInput potential
[out]vnewUpdated potential
[in]iterCurrent SCF iteration

Definition at line 1276 of file mix.F90.

◆ dbroyden_extrapolation()

subroutine mix_oct_m::dbroyden_extrapolation ( type(mix_t), intent(inout)  this,
type(namespace_t), intent(in)  namespace,
real(real64), intent(in)  coeff,
integer, intent(in)  d1,
integer, intent(in)  d2,
real(real64), dimension(:, :), intent(in)  vin,
real(real64), dimension(:, :), intent(out)  vnew,
integer, intent(in)  iter_used,
real(real64), dimension(:, :), intent(in)  f,
real(real64), dimension(:, :, :), intent(in)  df,
real(real64), dimension(:, :, :), intent(in)  dv 
)
private

Broyden mixing implementation.

Parameters
[in,out]thisMixing settings and derivatives
[in]coeffMixing coefficient
[in]d2Dimensions of the potential and residual
[in]iter_usedNumber of iterations to use in the mixing
[in]vinResidual and potential for current iteration
[in]dvResidual and potential from size(df, 4) iterations
[out]vnewUpdated potential

Definition at line 1323 of file mix.F90.

◆ dbroyden_extrapolation_aux()

subroutine mix_oct_m::dbroyden_extrapolation_aux ( type(mix_t), intent(inout)  this,
type(namespace_t), intent(in)  namespace,
integer, intent(in)  ii,
real(real64), intent(in)  coeff,
integer, intent(in)  iter_used,
real(real64), dimension(:, :), intent(in), optional  dbeta,
real(real64), dimension(:), intent(in), optional  dwork,
complex(real64), dimension(:, :), intent(in), optional  zbeta,
complex(real64), dimension(:), intent(in), optional  zwork 
)
private

Broyden mixing for auxilliary fields.

Auxilliary fields "follow" the direction of the mixing, but do not influence the direction. Auxilliary fields are passed as an attribute of a mix_t instance

Note
The design of this routine would be more consistent if mf (better-named aux_field) was passed in, then a wrapper method performs the iteration over the N aux fields, like:
do i = 1, this%nauxmixfield
aux_field => this%auxmixfield(i)%
if (field%func_type == type_float) then
call dbroyden_extrapolation_aux(this, namespace, aux_field, dbeta=beta, dwork=dwork)
else
call zbroyden_extrapolation_aux(this, namespace, aux_field, zbeta=beta, zwork=dwork)
end if
end do
rather than passing `ii` in.

Definition at line 1426 of file mix.F90.

◆ dmixing_diis()

subroutine mix_oct_m::dmixing_diis ( type(mix_t), intent(inout)  this,
real(real64), dimension(:, :), intent(in), contiguous  vin,
real(real64), dimension(:, :), intent(in), contiguous  residual,
real(real64), dimension(:, :), intent(out)  vnew,
integer, intent(in)  iter 
)
private

Direct inversion in the iterative subspace.

For implementation details, please refer to:

  • Pulay. [Convergence acceleration of iterative sequences the case of scf iteration](https:
  • G. Kresse, and J. Furthmueller. [Efficient iterative schemes for ab initio total-energy calculations using a plane-wave basis set](https:

Definition at line 1505 of file mix.F90.

◆ dmixing_build_matrix()

subroutine mix_oct_m::dmixing_build_matrix ( type(mix_t), intent(in)  this,
real(real64), dimension(:,:,:), intent(in)  df,
integer, intent(in)  size,
integer, intent(in)  d2,
real(real64), intent(in)  ww,
real(real64), dimension(:,:), intent(out), contiguous  beta 
)
private

Definition at line 1587 of file mix.F90.

◆ dmix_dotp()

real(real64) function mix_oct_m::dmix_dotp ( type(mix_t), intent(in)  this,
real(real64), dimension(:), intent(in)  xx,
real(real64), dimension(:), intent(in)  yy,
logical, intent(in), optional  reduce 
)
private

Definition at line 1641 of file mix.F90.

◆ dmixfield_set_vin()

subroutine mix_oct_m::dmixfield_set_vin ( type(mixfield_t), intent(inout)  mixfield,
real(real64), dimension(:, :), intent(in), contiguous  vin 
)
private

Definition at line 1678 of file mix.F90.

◆ dmixfield_set_vout()

subroutine mix_oct_m::dmixfield_set_vout ( type(mixfield_t), intent(inout)  mixfield,
real(real64), dimension(:, :), intent(in), contiguous  vout 
)
private

Definition at line 1690 of file mix.F90.

◆ dmixfield_get_vnew()

subroutine mix_oct_m::dmixfield_get_vnew ( type(mixfield_t), intent(in)  mixfield,
real(real64), dimension(:,:), intent(inout), contiguous  vnew 
)
private

Definition at line 1702 of file mix.F90.

◆ zmixing()

subroutine, public mix_oct_m::zmixing ( type(namespace_t), intent(in)  namespace,
type(mix_t), intent(inout), target  smix,
complex(real64), dimension(:, :), intent(in), contiguous  vin,
complex(real64), dimension(:, :), intent(in), contiguous  vout,
complex(real64), dimension(:, :), intent(out), contiguous  vnew 
)

Mix the input and output potentials (or densities) according to some scheme.

Whether \(v\) is the potential or the density is determined at run time.

Parameters
[in,out]smixMixing settings and derivatives
[in]vinInput potential
[in]voutOutput potential
[out]vnewUpdated potential

Definition at line 1793 of file mix.F90.

◆ zmixing_linear()

subroutine mix_oct_m::zmixing_linear ( integer, intent(in)  d1,
integer, intent(in)  d2,
real(real64), intent(in)  coeff,
complex(real64), dimension(:, :), intent(in)  vin,
complex(real64), dimension(:, :), intent(in)  residual,
complex(real64), dimension(:, :), intent(out)  vnew 
)
private

Linear mixing of the input and output potentials.

Linear mixing is defined as:

\[ V_{\text{new}} = (1 - a) V_{\text{in}} + a V_{\text{out}}. \]

Using the definition of the residual, f$ V_{\text{out}} = V_{\text{in}} + R

Parameters
[in]d2Limits for each dimension of the potential
[in]coeffMixing coefficient
[in]vinInput potential
[out]vnewUpdated potential

Definition at line 1847 of file mix.F90.

◆ zcompute_residual()

subroutine mix_oct_m::zcompute_residual ( type(mix_t), intent(in)  smix,
integer, dimension(2), intent(in)  dims,
complex(real64), dimension(:, :), intent(in)  vin,
complex(real64), dimension(:, :), intent(in)  vout,
complex(real64), dimension(:, :), intent(out), contiguous  residual 
)
private

Compute the residual of the potential (or density) for SCF mixing.

In the absence of any preconditioning, the residual is defined as \( V_\text{out} - V_\text{in}\). If preconditioning is used, please refer to the specific scheme for the definition.

Note, d1, d2 are passed in because this routine should be compatible with both the potential and and the auxiliary fields (which have differing sizes).

Parameters
[in]smixMixing settings and derivatives
[in]dimsDimensions of the terms to mix
[in]vinInput potential
[in]voutOutput potential

Definition at line 1879 of file mix.F90.

◆ zmixing_broyden()

subroutine mix_oct_m::zmixing_broyden ( type(namespace_t), intent(in)  namespace,
type(mix_t), intent(inout)  smix,
complex(real64), dimension(:, :), intent(in), contiguous  f,
complex(real64), dimension(:, :), intent(in), contiguous  vin,
complex(real64), dimension(:, :), intent(out), contiguous  vnew,
integer, intent(in)  iter 
)
private

Top-level routine for Broyden mixing.

Performs updating of smixmixfield attributes. Whether \(v\) is the potential or the density is determined at run time.

Parameters
[in,out]smixMixing settings and derivatives
[in]fResidual
[in]vinInput potential
[out]vnewUpdated potential
[in]iterCurrent SCF iteration

Definition at line 1929 of file mix.F90.

◆ zbroyden_extrapolation()

subroutine mix_oct_m::zbroyden_extrapolation ( type(mix_t), intent(inout)  this,
type(namespace_t), intent(in)  namespace,
real(real64), intent(in)  coeff,
integer, intent(in)  d1,
integer, intent(in)  d2,
complex(real64), dimension(:, :), intent(in)  vin,
complex(real64), dimension(:, :), intent(out)  vnew,
integer, intent(in)  iter_used,
complex(real64), dimension(:, :), intent(in)  f,
complex(real64), dimension(:, :, :), intent(in)  df,
complex(real64), dimension(:, :, :), intent(in)  dv 
)
private

Broyden mixing implementation.

Parameters
[in,out]thisMixing settings and derivatives
[in]coeffMixing coefficient
[in]d2Dimensions of the potential and residual
[in]iter_usedNumber of iterations to use in the mixing
[in]vinResidual and potential for current iteration
[in]dvResidual and potential from size(df, 4) iterations
[out]vnewUpdated potential

Definition at line 1976 of file mix.F90.

◆ zbroyden_extrapolation_aux()

subroutine mix_oct_m::zbroyden_extrapolation_aux ( type(mix_t), intent(inout)  this,
type(namespace_t), intent(in)  namespace,
integer, intent(in)  ii,
real(real64), intent(in)  coeff,
integer, intent(in)  iter_used,
real(real64), dimension(:, :), intent(in), optional  dbeta,
real(real64), dimension(:), intent(in), optional  dwork,
complex(real64), dimension(:, :), intent(in), optional  zbeta,
complex(real64), dimension(:), intent(in), optional  zwork 
)
private

Broyden mixing for auxilliary fields.

Auxilliary fields "follow" the direction of the mixing, but do not influence the direction. Auxilliary fields are passed as an attribute of a mix_t instance

Note
The design of this routine would be more consistent if mf (better-named aux_field) was passed in, then a wrapper method performs the iteration over the N aux fields, like:
do i = 1, this%nauxmixfield
aux_field => this%auxmixfield(i)%
if (field%func_type == type_float) then
call dbroyden_extrapolation_aux(this, namespace, aux_field, dbeta=beta, dwork=dwork)
else
call zbroyden_extrapolation_aux(this, namespace, aux_field, zbeta=beta, zwork=dwork)
end if
end do
rather than passing `ii` in.

Definition at line 2079 of file mix.F90.

◆ zmixing_diis()

subroutine mix_oct_m::zmixing_diis ( type(mix_t), intent(inout)  this,
complex(real64), dimension(:, :), intent(in), contiguous  vin,
complex(real64), dimension(:, :), intent(in), contiguous  residual,
complex(real64), dimension(:, :), intent(out)  vnew,
integer, intent(in)  iter 
)
private

Direct inversion in the iterative subspace.

For implementation details, please refer to:

  • Pulay. [Convergence acceleration of iterative sequences the case of scf iteration](https:
  • G. Kresse, and J. Furthmueller. [Efficient iterative schemes for ab initio total-energy calculations using a plane-wave basis set](https:

Definition at line 2158 of file mix.F90.

◆ zmixing_build_matrix()

subroutine mix_oct_m::zmixing_build_matrix ( type(mix_t), intent(in)  this,
complex(real64), dimension(:,:,:), intent(in)  df,
integer, intent(in)  size,
integer, intent(in)  d2,
real(real64), intent(in)  ww,
complex(real64), dimension(:,:), intent(out), contiguous  beta 
)
private

Definition at line 2240 of file mix.F90.

◆ zmix_dotp()

complex(real64) function mix_oct_m::zmix_dotp ( type(mix_t), intent(in)  this,
complex(real64), dimension(:), intent(in)  xx,
complex(real64), dimension(:), intent(in)  yy,
logical, intent(in), optional  reduce 
)
private

Definition at line 2294 of file mix.F90.

◆ zmixfield_set_vin()

subroutine mix_oct_m::zmixfield_set_vin ( type(mixfield_t), intent(inout)  mixfield,
complex(real64), dimension(:, :), intent(in), contiguous  vin 
)
private

Definition at line 2331 of file mix.F90.

◆ zmixfield_set_vout()

subroutine mix_oct_m::zmixfield_set_vout ( type(mixfield_t), intent(inout)  mixfield,
complex(real64), dimension(:, :), intent(in), contiguous  vout 
)
private

Definition at line 2343 of file mix.F90.

◆ zmixfield_get_vnew()

subroutine mix_oct_m::zmixfield_get_vnew ( type(mixfield_t), intent(in)  mixfield,
complex(real64), dimension(:,:), intent(inout), contiguous  vnew 
)
private

Definition at line 2355 of file mix.F90.

Variable Documentation

◆ max_auxmixfield

integer, parameter mix_oct_m::max_auxmixfield = 5
private

Definition at line 204 of file mix.F90.