solnScatPot
index

potential term for solution-phase scattering experiments, including 
SAXS and SANS
 
 One normally creates SolnScatPot objects using 
solnXRayPotTools.create_XRayScatPot or sansPotTools.create_SANSPot.
 
This term is ensembleSimulation.EnsembleSimulation - aware.
 
See C.D. Schwieters and G.M. Clore, ``Using Small Angle Solution
          Scattering Data in Xplor-NIH Structure Calculations ,''
          Progr. NMR Spectroscopy, accepted (2014).
 
constructor:
 
  SolnScatPot(instanceName ,
              sel          ,
              qValues      ,
              groupMap     ,
              formFactors  ,
              iFormFactors ,
              solventParams,
              rho0           ,         
              radiusScale  ,
              volumeScale  ,
              globDefs     ,
              experiment   ,
              qValuesExpt  ,
              simulation   )
 
  where:
    instanceName  user-specified identifier.
    sel           atomSel.AtomSel which specifies atoms to use in 
                  the scattering calculation (frequently all heavy atoms).
    qValues       a sequence of values of scattering vector at which
                  to perform refinement. q=4 PI sin(theta) / lambda,
                  where 2 theta is the scattering angle, and lambda is
                  the wavelength of the incident beam. These should usually
                  be specified at uniform values of q.
    groupMap      a sequence of group names to use in the scattering calc.
    formFactors   a sequence of maps of group name to effective form factor,
                  one for each value of scattering vector- in the same order
    iFormFactors  a sequence of maps of atom index to effective form factor,
                  one for each value of scattering vector- in the same order.
                  Generally, individual atoms need not be specified-
                  but are rather specified by group name. This is for 
                  exceptions, and will override any group entry.
    solventParams a dictionary of 2-membered tuples- one for each atom group
                  name. Each tuple corresponds to the radius and volume for 
                  that atom group.
    rho0          excluded solvent electron density.
    radiusScale   scale factor used in the calculation of the displaced
                  solvent scattering contribution. See below.
 
    volumeScale   scale factor used in the calculation of the displaced
                  solvent scattering contribution. See below.
    globDefs      a sequence of glob definitions. Each definition is a 
                  sequence of tuples (w,index), where w is a weight and index
                  is an atom index.
    qValuesExpt   The q-values corresponding to the experimental measurements.
                  The fit of calculated curve to experiment is done for
                  each value of q here, with the calculated curve interpolated 
                  using a cubic spline, hence it is important that none of
                  the qValuesExpt fall outside the range given in qValues.
    experiment    solvent-corrected experimental scattering intensity values-
                  one for each value of q in qValuesExpt.
    simulation    optional simulation.Simulation object.
 
members:
 The following parameters can be set [defaults in square brackets]
 
  scale            - scale factor (force constant) [1]
 
  weights   - sequence of factors to weight calculated scattering intensity
              to its corresponding experimental value. Used for rms, energy 
              calculation. See below. Default is a value of 1 for all.
 
  calcType  - specifies how I(q) is to be calculated. Valid values are
                "n2"  - use the exact, but slow Debye formula.
                "n"   - numerically integrate over surface of constant
                        q, using theta, phi values determined by angleType
                 "uniform" - a faster version of the n algorithm which
                             can be used if q values are equally spaced.
                 "multipole" - Use the spherical harmonic expansion of
                                 Svergun, Barberato and Koch, 
                                 J. Appl. Cryt. 28, 768 (1995)
                               Note: at this time it is not possible
                               to refine using the calcType, as the 
                               derivatives are not yet coded.
 
 
  numAngles - for calcTypes n and uniform, the number of (theta,phi)
              pairs to use in the evaluation of I(q).
 
  qValues   - the values of scattering vector amplitude q.
 
  expt()    - the observed values of I(q).
 
  formFactors  - specified the formFactor map.
 
  iFormFactors - specify the atom-indexed formFactor map.
 
  angleType - specifies how (theta,phi) angle pairs are chosen over
              the surface of a sphere. Choices are:
                "random"    - choose randomly- not recommended.
                "fibonacci" - method based on Fibonacci numbers.
                                D. Svergun, Acta Cryst. A50, 391 (1994).
                "spiral"    - spiral algorithm from 
                                E.B. Saff and A.B.J. Kuijlaars, 
                                Distributing Many Points on a Sphere,
                                The Mathematical Intelligencer, 19(1), 
                                Winter (1997).
                              This is the preferred algorithm.
 
  lMax      - for calcType multipole, the maximum value of l to use in
              the spherical harmonic expansion.
 
  cmpType   - should be "plain" or "log." This controls whether the
              logarithm is taken before the energy is evaluated. See below.
              [Default: plain]
 
  globs     - sequence of glob definitions - see description above.
 
  normalizeIndex - if greater than or equal to zero, specifies which
                   value of q at which to consider Icalc and Iobs to be
                   equal. If the value is -2, the normalization is taken
                   as the average of the Icalc or Iobs values weighted by
                   weights. If the value of normalizeIndex lies outside 
                   the above values, no normalization is performed. If the
                   value -3 is specified, the normalization will instead 
                   minimize chi^2. 
                   [default: -3]
 
  radiusScale - scale factor used in the calculation of the displaced
                solvent scattering contribution. See below.
 
  volumeScale - scale factor used in the calculation of the displaced
                solvent scattering contribution. See below.
 
  radiusType  - "radius" or "volume" - specifies the method used for 
                computing atomic radii. See below.
 
  rho0        - excluded solvent electron density.
 
  calcBoundary         - boolean - whether to calculate boundary layer 
                         scattering contribution [default: False]
 
  rhob                 - boundary layer electron density [default: 0.03]
 
  boundaryThickness    - the thickness of the boundary layer [default: 2.88].
 
  solventRadius        - radius used in the surface computation 
                         [default: 1.44]. 
 
  bg                   - isotropic background value to be subtracted from
                         observed scattering values.
 
  formScale   - overall scale factor for vacuum atomic form factors.
 
  useGlobs  - whether or not to use globs in the scattering calculation.
 
  verbose   - if true, produce more detailed output [False]
 
  useSimEnsWeight - whether to use the ensemble wieghts set with setEnsWeights
                    or to use those of the underlying EnsembleSimulation.
  ensWeights      - a sequence of ensemble weights to use when calculating 
                    the ensemble- averaged value of each scattering curve. By 
                    default, the weights of the underlying EnsembleSimulation 
                    are used. If these are overridden by calling the 
                    setEnsWeights method, and useSimEnsWeight is set to False. 
 
 
 the above quantities may be retrieved using the member function form
 quantity(), while they are set using the form setQuantity(value).
 
 
  QValues - a plain member containing values of (x,y,z) on a r=1
            sphere used for numerical integration of I(q).
 
 
methods:
  
  calcEnergy()                 - calc energy, returns the energy value.
  calcEnergyAndDerivs(derivs)  - calc energy, derivs, returns the energy value.
 
  calcGlobCorrect(val) - calculate the globic correction factor. The argument
                         val can be calcType (e.g. 'n2'), in which case 
                         the globic correction factor will be calculated by
                         using the specified calcType. val can alternatively
                         be a raw array of I(q) values which are used to
                         generate g(q). val defaults to the current calcType.
 
  globCorrect() - return the current globic correction factor.
 
  getF()        - return a matrix of scattering amplitudes at each angle
                  and q value.
 
  rms()         - return the rms of calcedShift - effShift
 
  info()       - string with current info about the state of this instance
  simulation() - return the associated simulation.Simulation.
 
  calcd()      - return the calculated values of I(q) for each value of 
                 qValues.
  splined()    - return the calculated values of I(q) for each value of 
                 qValuesExpt.
 
  exptScale()  Iobs(qn), for n=normalizeIndex, if normalizeIndex>-1, 
               1 otherwise.
  calcScale()  Icalc(qn), for n=normalizeIndex, if normalizeIndex>-1, 
               1 otherwise.
 
  I_contrib(memberIndex) - the scatter intensity from the specified ensemble
                           member structure, or the local member if the
                           index is omitted.
  boundaryI()            - scattering intensity contribution from the boundary
                           layer.
  selection()            - return the atom selection specified in the 
                           constructor.
  radii()                - return array of atomic radii
  aveRadius()            - average atomic radius.
  aveVolume()            - average atomic volume.
 
  boundaryVol()          - volume of the boundary layer.
  
  calcBoundaryPoints()   - calculate new points used to compute the 
                           boundary layer scattering contribution.
 
  ensWeight(index)       - return the ensemble weight associated with the
                           specified member.
 
  addRigidRegion(atomSel) - indicate that the specified atoms move as
                            a rigid body so that a very fast update
                            algorithm can be used. Use of this method
                            does not actually constrain the relative
                            atom positions to be constrained. For
                            that, the atoms must be grouped in the
                            appropriate ivm.IVM objects.
 
  solnScat()              - return a SolnScat object. This object is
                            used to calculate the raw scattering curve
                            from atomic coordinates. It has a few
                            useful methods to access properties not
                            directly accessible from the SolnScatPot
                            class:
 
      numRigidRegions()       - return the number of rigid regions.
      rigidRegion( n )        - return the definition of rigid region
                                n, indexed from 0.
 
        Accessor members:
      innerBoundaryFilename   - if set, write out the inner molecular surface
                                used for boundary contribution calculation.
      outerBoundaryFilename   - if set, write out the outer molecular surface
                                used for boundary contribution calculation.
 
Energy is calculated as
 
  scale * rms^2
 
where
 
   rms^2 = sum_i w_i (Icalc_i - Iobs_i)^2 / (Nq),
 
in which Nq is the total number of values of q at which Iobs is
specified, and w_i is the corresponding weight. Typically, we choose
w_i = (sigma_i)^(-2), where sigma_i is the error in I(q) at q=q_i.
 
Given a set of atomic coordinates the scattering intensity is
calculated from the following expression: 
 
    I(q) = <| A(Q) |^2>_q
 
where <>_q denotes average over solid angle at constant amplitude of
the scattering vector.
 
The scattering amplitude A(Q) consists of a sum of contributions from
each atom (usually only heavy atoms are included)
 
  A(Q) = sum_j f_j(q) exp(i dot(Q,r_j)) + 
         sum_j2 rhob(r_j2) exp(i dot(Q,r_j2))
 
where f_j(q) is the total effective scattering amplitude of atom j
(including excluded solvent effects), and r_j is atom j's position. The
second contribution is due to boundary-layer scattering, and the sum
over j2 is, in principle an integral over the region in which rhob,
the boundary layer density is nonzero.
 
For rhob=0, the integral over solid angle can be carried out
analytically, resulting in the Debye formula calcType=n2:
 
  I(q) = sum_j1 sum_j2 f_j1(q) f_j2(q) sinc(q |r_j1-r_j2|)
 
where sinc(x) = sin(x) / x.
 
This expression is quite expensive for biological systems due to the
double sum.
 
One alternative is to express A(Q) as an expansion in terms of
spherical harmonic functions (calcType=multipole), and approximate
I(q) by truncating the expansion. This is works well (is efficient)
for small scattering angles - like those used in most SAXS
experiments. For larger angles, the expansion fails spectacularly.
 
One can also evaluate the integral over solid angle by simple grid
integration (calcTypes=n, uniform). This works well for all moderate
values of q.
 
Excluded Solvent Treatment:
 
  The total atomic effective scattering amplitude is
 
   f_j(q) = fv_j(q) - g_j(q),
 
where fv_j(q) is the vacuum atomic form factor, and g_j is an excluded
solvent contribution
 
  g_j(q) = rho0 * V_j * exp(-s**2 * PI * V_j**(2./3)) *
           volumeScale * exp(-q**2 * PI * (4PI/3)^{1/3} * (r0^2 - rm^2)),
 
which depends on atomic volume V_j and scaled scattering vector
amplitude s=q/(2PI). rm is aveRadius() and r0 is a scaled radius
 
   r0= radiusScale() * rm
 
If radiusType=="volume" (the default), rm actually is computed as the
radius corresoponding to aveVolume() [i.e. (3 Vm / (4 PI))^(1/3) ]
 
Boundary Layer Treatment:
 
The boundary layer contribution is computed based on a molecular
surface generated using surfD. The molecular surface is generated
based on atomic radii, boundaryThickness and solventRadius(). An outer
surface is generated by rolling a ball of radius solventRadius over
the surface generated of spheres of radius
radii+boundaryThickness. Then, an inner surface is generated by, for
each surface triangle, dropping a line segment of length
boundaryThickness in the direction opposite the surface normal. Each
surface triangle thus generates a voxel. The contribution of each
voxel to the scattering amplitude is represented by a sphere sharing
the center and volume of the voxel.
 
Ensemble Calculations
 
When calculated in the context of an ensembleSimulation.EnsembleSimulation
the scattering intensity becomes:
 
   I(q) = sum_n W_n I_n(q)
 
where the sum is over all members of an ensemble, w_n is the weight
for ensemble member n, and I_n(q) is the scattering intensity from
ensemble member n, calculated as above.
 
One should note that the approximate methods of calculating I(q)
introduce a rotational dependence: the grid on which the scattering
amplitude is calculated is fixed in space. This rotational dependence
implies that the potential energy term is not constant upon overall
rotation. The amount of the rotational dependence decreases with
increased numAngles. A future enhancement may have the grid rotate
along with the rest of the system.
 
 
 
# This file was automatically generated by SWIG (http://www.swig.org).
# Version 4.0.2
#
# Do not make changes to this file unless you know what you are doing--modify
# the SWIG interface file instead.

 
Classes
       
builtins.object
CDSList_DComplex
CDSList_EigenPair_DComplex
CDSList_EigenPair_double
CDSList_FloatPair
CDSList_Voxel
CDSMatrix_complex
CDSVector_CDSVector_Complex
CDSVector_complex
EigenPair_DComplex
EigenPair_double
EnsemblePot
SolnScatPot_LetterClass
InterpSphere_CDSVector_complex
Modified
ModifiedBase
SolnScat
SolnScatPot
SolnScatRigid
SolnScat_FracElec_SelIndex
VarEnsWeights
Voxel
rc_EnsemblePot

 
class CDSList_DComplex(builtins.object)
    CDSList_DComplex(*args)
 

 
  Methods defined here:
__delitem__(self, *args, **kwargs) -> 'void'
__getitem__(self, *args) -> 'CDSList< DComplex >'
__getslice__(self, *args, **kwargs) -> 'CDSList< DComplex >'
__init__(self, *args)
Initialize self.  See help(type(self)) for accurate signature.
__len__(self, *args, **kwargs) -> 'int'
__repr__ = _swig_repr(self)
__setitem__(self, *args, **kwargs) -> 'void'
append(self, *args, **kwargs) -> 'void'
help(self, *args, **kwargs) -> 'String'
remove(self, *args, **kwargs) -> 'void'
removeAll(self, *args, **kwargs) -> 'void'

Static methods defined here:
__swig_destroy__ = delete_CDSList_DComplex(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag

 
class CDSList_EigenPair_DComplex(builtins.object)
    CDSList_EigenPair_DComplex(*args)
 

 
  Methods defined here:
__delitem__(self, *args, **kwargs) -> 'void'
__getitem__(self, *args) -> 'CDSList< EigenPair< DComplex > >'
__getslice__(self, *args, **kwargs) -> 'CDSList< EigenPair< DComplex > >'
__init__(self, *args)
Initialize self.  See help(type(self)) for accurate signature.
__len__(self, *args, **kwargs) -> 'int'
__repr__ = _swig_repr(self)
__setitem__(self, *args, **kwargs) -> 'void'
append(self, *args, **kwargs) -> 'void'
help(self, *args, **kwargs) -> 'String'
remove(self, *args, **kwargs) -> 'void'
removeAll(self, *args, **kwargs) -> 'void'

Static methods defined here:
__swig_destroy__ = delete_CDSList_EigenPair_DComplex(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag

 
class CDSList_EigenPair_double(builtins.object)
    CDSList_EigenPair_double(*args)
 

 
  Methods defined here:
__delitem__(self, *args, **kwargs) -> 'void'
__getitem__(self, *args) -> 'CDSList< EigenPair< double > >'
__getslice__(self, *args, **kwargs) -> 'CDSList< EigenPair< double > >'
__init__(self, *args)
Initialize self.  See help(type(self)) for accurate signature.
__len__(self, *args, **kwargs) -> 'int'
__repr__ = _swig_repr(self)
__setitem__(self, *args, **kwargs) -> 'void'
append(self, *args, **kwargs) -> 'void'
help(self, *args, **kwargs) -> 'String'
remove(self, *args, **kwargs) -> 'void'
removeAll(self, *args, **kwargs) -> 'void'

Static methods defined here:
__swig_destroy__ = delete_CDSList_EigenPair_double(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag

 
class CDSList_FloatPair(builtins.object)
    CDSList_FloatPair(*args)
 

 
  Methods defined here:
__delitem__(self, *args, **kwargs) -> 'void'
__getitem__(self, *args) -> 'CDSList< Pair< double,double > >'
__getslice__(self, *args, **kwargs) -> 'CDSList< Pair< double,double > >'
__init__(self, *args)
Initialize self.  See help(type(self)) for accurate signature.
__len__(self, *args, **kwargs) -> 'int'
__repr__ = _swig_repr(self)
__setitem__(self, *args, **kwargs) -> 'void'
append(self, *args, **kwargs) -> 'void'
help(self, *args, **kwargs) -> 'String'
remove(self, *args, **kwargs) -> 'void'
removeAll(self, *args, **kwargs) -> 'void'

Static methods defined here:
__swig_destroy__ = delete_CDSList_FloatPair(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag

 
class CDSList_Voxel(builtins.object)
    CDSList_Voxel(*args)
 

 
  Methods defined here:
__delitem__(self, *args, **kwargs) -> 'void'
__getitem__(self, *args) -> 'CDSList< Voxel >'
__getslice__(self, *args, **kwargs) -> 'CDSList< Voxel >'
__init__(self, *args)
Initialize self.  See help(type(self)) for accurate signature.
__len__(self, *args, **kwargs) -> 'int'
__repr__ = _swig_repr(self)
__setitem__(self, *args, **kwargs) -> 'void'
append(self, *args, **kwargs) -> 'void'
help(self, *args, **kwargs) -> 'String'
remove(self, *args, **kwargs) -> 'void'
removeAll(self, *args, **kwargs) -> 'void'

Static methods defined here:
__swig_destroy__ = delete_CDSList_Voxel(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag

 
class CDSMatrix_complex(builtins.object)
    CDSMatrix_complex(*args)
 

 
  Methods defined here:
__add__(self, *args, **kwargs) -> 'CDSMatrix< CDS::Complex< double > >'
__getitem__(self, *args, **kwargs) -> 'CDS::Complex< double >'
__iadd__(self, *args, **kwargs) -> 'CDSMatrix< CDS::Complex< double > >'
__init__(self, *args)
Initialize self.  See help(type(self)) for accurate signature.
__mul__(self, *args) -> 'CDSMatrix< CDS::Complex< double > >'
__neg__(self, *args, **kwargs) -> 'CDSMatrix< CDS::Complex< double > >'
__repr__(self, *args, **kwargs) -> 'String'
Return repr(self).
__rmul__(self, *args, **kwargs) -> 'CDSMatrix< CDS::Complex< double > >'
__setitem__(self, *args, **kwargs) -> 'void'
__str__(self, *args, **kwargs) -> 'String'
Return str(self).
__sub__(self, *args, **kwargs) -> 'CDSMatrix< CDS::Complex< double > >'
add(self, *args, **kwargs) -> 'void'
cols(self, *args, **kwargs) -> 'int'
fromList(s, l)
get(self, *args, **kwargs) -> 'CDS::Complex< double >'
rawData(self, *args, **kwargs) -> 'CDSVector< CDS::Complex< double > >'
resize(self, *args, **kwargs) -> 'void'
rows(self, *args, **kwargs) -> 'int'
scale(self, *args, **kwargs) -> 'void'
set(self, *args) -> 'void'
setDiag(self, *args) -> 'void'

Static methods defined here:
__swig_destroy__ = delete_CDSMatrix_complex(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag

 
class CDSVector_CDSVector_Complex(builtins.object)
    CDSVector_CDSVector_Complex(*args)
 

 
  Methods defined here:
__add__(self, *args) -> 'CDSVector< CDSVector< CDS::complex > >'
__call__(self, *args) -> 'CDSVector< CDS::Complex< double >,0,CDS::DefaultAlloc,long,int > const &'
Call self as a function.
__getitem__(self, *args, **kwargs)
__getslice__(self, *args, **kwargs) -> 'CDSVector< CDSVector< CDS::complex > >'
__iadd__(self, *args) -> 'CDSVector< CDSVector< CDS::complex > >'
__imul__(self, *args, **kwargs) -> 'CDSVector< CDSVector< CDS::complex > >'
__init__(self, *args)
Initialize self.  See help(type(self)) for accurate signature.
__isub__(self, *args) -> 'CDSVector< CDSVector< CDS::complex > >'
__itruediv__(self, *args, **kwargs) -> 'CDSVector< CDSVector< CDS::complex > >'
__len__(self, *args, **kwargs) -> 'long'
__mul__(self, *args) -> 'CDSVector< CDSVector< CDS::complex > >'
__neg__(self, *args, **kwargs) -> 'CDSVector< CDSVector< CDS::complex > >'
__pow__(self, *args, **kwargs) -> 'CDSVector< CDSVector< CDS::complex > >'
__radd__(self, *args, **kwargs) -> 'CDSVector< CDSVector< CDS::complex > >'
__repr__ = _swig_repr(self)
__rmul__(self, *args, **kwargs) -> 'CDSVector< CDSVector< CDS::complex > >'
__rsub__(self, *args, **kwargs) -> 'CDSVector< CDSVector< CDS::complex > >'
__rtruediv__(self, *args, **kwargs) -> 'CDSVector< CDSVector< CDS::complex > >'
__setitem__(self, *args, **kwargs) -> 'void'
__str__(self, *args, **kwargs) -> 'String'
Return str(self).
__sub__(self, *args) -> 'CDSVector< CDSVector< CDS::complex > >'
__truediv__(self, *args) -> 'CDSVector< CDSVector< CDS::complex > >'
append(self, *args, **kwargs) -> 'void'
fromList(s, l)
get(self, *args, **kwargs) -> 'CDSVector< CDS::Complex< double >,0,CDS::DefaultAlloc,long,int >'
help(self, *args, **kwargs) -> 'String'
offset(self, *args, **kwargs) -> 'int'
scale(self, *args, **kwargs) -> 'void'

Static methods defined here:
__swig_destroy__ = delete_CDSVector_CDSVector_Complex(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag

 
class CDSVector_complex(builtins.object)
    CDSVector_complex(*args)
 

 
  Methods defined here:
__add__(self, *args) -> 'CDSVector< CDS::complex >'
__call__(self, *args) -> 'CDS::Complex< double > const &'
Call self as a function.
__getitem__(self, *args, **kwargs)
__getslice__(self, *args, **kwargs) -> 'CDSVector< CDS::complex >'
__iadd__(self, *args) -> 'CDSVector< CDS::complex >'
__imul__(self, *args, **kwargs) -> 'CDSVector< CDS::complex >'
__init__(self, *args)
Initialize self.  See help(type(self)) for accurate signature.
__isub__(self, *args) -> 'CDSVector< CDS::complex >'
__itruediv__(self, *args, **kwargs) -> 'CDSVector< CDS::complex >'
__len__(self, *args, **kwargs) -> 'long'
__mul__(self, *args) -> 'CDSVector< CDS::complex >'
__neg__(self, *args, **kwargs) -> 'CDSVector< CDS::complex >'
__pow__(self, *args, **kwargs) -> 'CDSVector< CDS::complex >'
__radd__(self, *args, **kwargs) -> 'CDSVector< CDS::complex >'
__repr__ = _swig_repr(self)
__rmul__(self, *args, **kwargs) -> 'CDSVector< CDS::complex >'
__rsub__(self, *args, **kwargs) -> 'CDSVector< CDS::complex >'
__rtruediv__(self, *args, **kwargs) -> 'CDSVector< CDS::complex >'
__setitem__(self, *args, **kwargs) -> 'void'
__str__(self, *args, **kwargs) -> 'String'
Return str(self).
__sub__(self, *args) -> 'CDSVector< CDS::complex >'
__truediv__(self, *args) -> 'CDSVector< CDS::complex >'
append(self, *args, **kwargs) -> 'void'
fromList(s, l)
get(self, *args, **kwargs) -> 'CDS::Complex< double >'
help(self, *args, **kwargs) -> 'String'
offset(self, *args, **kwargs) -> 'int'
scale(self, *args, **kwargs) -> 'void'

Static methods defined here:
__swig_destroy__ = delete_CDSVector_complex(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag

 
class EigenPair_DComplex(builtins.object)
    EigenPair_DComplex(*args, **kwargs)
 

 
  Methods defined here:
__init__(self, *args, **kwargs)
Initialize self.  See help(type(self)) for accurate signature.
__repr__ = _swig_repr(self)
value(self, *args, **kwargs) -> 'DComplex'
vector(self, *args, **kwargs) -> 'CDSList< DComplex >'

Static methods defined here:
__swig_destroy__ = delete_EigenPair_DComplex(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag
value_

 
vector_

 

 
class EigenPair_double(builtins.object)
    EigenPair_double(*args, **kwargs)
 

 
  Methods defined here:
__init__(self, *args, **kwargs)
Initialize self.  See help(type(self)) for accurate signature.
__repr__ = _swig_repr(self)
value(self, *args, **kwargs) -> 'double'
vector(self, *args, **kwargs) -> 'CDSList< double >'

Static methods defined here:
__swig_destroy__ = delete_EigenPair_double(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag
value_

 
vector_

 

 
class EnsemblePot(builtins.object)
    EnsemblePot(*args, **kwargs)
 

 
  Methods defined here:
__init__(self, *args, **kwargs)
Initialize self.  See help(type(self)) for accurate signature.
__repr__ = _swig_repr(self)
addEnsWeights(self, *args, **kwargs) -> 'void'
calcEnergy(self, *args, **kwargs) -> 'float_type'
calcEnergyAndDerivs(self, *args, **kwargs) -> 'float_type'
calcWDerivs(self, *args, **kwargs) -> 'bool const'
clearEnsWeights(self, *args, **kwargs) -> 'void'
energyMaybeDerivs0(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivs1(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivs2(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivs3(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivs4(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivsPost(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivsPre(self, *args, **kwargs) -> 'float_type'
ensWeight(self, *args, **kwargs) -> 'float_type'
ensWeights(self, *args, **kwargs) -> 'CDSList< float_type >'
ensWeightsInfo(self, *args, **kwargs) -> 'String'
getEnsWeights(self, *args, **kwargs) -> 'CDSList< VarEnsWeights > &'
setCalcWDerivs(self, *args, **kwargs) -> 'void'
setEnsWeights(self, *args, **kwargs) -> 'void'
setUseSimEnsWeights(self, *args, **kwargs) -> 'void'
simulation(self, *args) -> 'EnsembleSimulationBase const *'
updateEnsWeights(self, *args, **kwargs) -> 'void'
useSimEnsWeights(self, *args, **kwargs) -> 'bool const'

Static methods defined here:
__swig_destroy__ = delete_EnsemblePot(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag

 
class InterpSphere_CDSVector_complex(builtins.object)
    InterpSphere_CDSVector_complex(*args, **kwargs)
 

 
  Methods defined here:
__call__(self, *args, **kwargs) -> 'CDSVector< CDS::Complex< double >,0,CDS::DefaultAlloc,long,int >'
Call self as a function.
__init__(self, *args, **kwargs)
Initialize self.  See help(type(self)) for accurate signature.
__repr__ = _swig_repr(self)
derivs(self, *args, **kwargs) -> 'XplorNIH::InterpSphere< CDSVector< CDS::Complex< double > > >::GradType'
valAndDerivs(self, *args, **kwargs) -> 'void'

Static methods defined here:
__swig_destroy__ = delete_InterpSphere_CDSVector_complex(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag

 
class Modified(builtins.object)
    Modified(*args, **kwargs)
 

 
  Methods defined here:
__call__(self, *args, **kwargs) -> 'int'
Call self as a function.
__init__(self, *args, **kwargs)
Initialize self.  See help(type(self)) for accurate signature.
__repr__ = _swig_repr(self)
asString(self, *args, **kwargs) -> 'String'
clear(self, *args, **kwargs) -> 'void'
set(self, *args, **kwargs) -> 'void'
update(self, *args, **kwargs) -> 'void'
value(self, *args, **kwargs) -> 'int'

Static methods defined here:
__swig_destroy__ = delete_Modified(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag

Data and other attributes defined here:
MOD_SELF = 1
MOD_SIMULATION = 2

 
class ModifiedBase(builtins.object)
    ModifiedBase(*args, **kwargs)
 

 
  Methods defined here:
__init__(self, *args, **kwargs)
Initialize self.  See help(type(self)) for accurate signature.
__repr__ = _swig_repr(self)
registerTo(self, *args, **kwargs) -> 'void'
unRegister(self, *args, **kwargs) -> 'void'
updateValues(self, *args, **kwargs) -> 'void'

Static methods defined here:
__swig_destroy__ = delete_ModifiedBase(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
modified

 
registeredSimulations

 
thisown

 
The membership flag

 
class SolnScat(builtins.object)
    SolnScat(*args, **kwargs)
 

 
  Methods defined here:
__init__(self, *args, **kwargs)
Initialize self.  See help(type(self)) for accurate signature.
__repr__ = _swig_repr(self)
addRigidRegion(self, *args, **kwargs) -> 'void'
angleType(self, *args, **kwargs) -> 'SolnScat::SolnScatAngleType'
aveRadius(self, *args, **kwargs) -> 'float_type'
aveVolume(self, *args, **kwargs) -> 'float_type'
boundaryI(self, *args, **kwargs) -> 'CDSVector< float_type >'
boundaryThickness(self, *args, **kwargs) -> 'float_type const'
boundaryVol(self, *args, **kwargs) -> 'float_type'
calcBoundary(self, *args, **kwargs) -> 'bool const'
calcBoundaryPoints(self, *args, **kwargs) -> 'void'
calcDerivs(self, *args, **kwargs) -> 'CDSVector< Vec3 >'
calcFGlob(self, *args, **kwargs) -> 'void'
calcI(self, *args, **kwargs) -> 'CDSVector< float_type >'
calcIFromPos(self, *args, **kwargs) -> 'CDSVector< float_type >'
calcI_N(self, *args, **kwargs) -> 'CDSVector< float_type >'
calcI_N2(self, *args, **kwargs) -> 'CDSVector< float_type >'
calcI_multipole(self, *args, **kwargs) -> 'CDSVector< float_type >'
calcI_uniform(self, *args, **kwargs) -> 'CDSVector< float_type >'
calcPosGlob(self, *args, **kwargs) -> 'void'
calcType(self, *args, **kwargs) -> 'SolnScat::SolnScatCalcType const'
calc_g(self, *args, **kwargs) -> 'float_type'
convertFormFactors(self, *args, **kwargs) -> 'void'
formFactors(self, *args, **kwargs) -> 'CDSVector< SolnScat::FormFactorMap > const &'
formScale(self, *args, **kwargs) -> 'float_type'
freeSurfVoxels(self, *args, **kwargs) -> 'bool const'
getF(self, *args, **kwargs) -> 'CDSMatrix< CDS::complex > const'
globs(self, *args, **kwargs) -> 'CDSList< GlobDef > const &'
iFormFactors(self, *args, **kwargs) -> 'CDSVector< SolnScat::FormFactorIMap > const &'
innerBoundaryFilename(self, *args, **kwargs) -> 'String const'
lMax(self, *args, **kwargs) -> 'int const'
numAngles(self, *args, **kwargs) -> 'int'
numRigidRegions(self, *args, **kwargs) -> 'int'
outerBoundaryFilename(self, *args, **kwargs) -> 'String const'
qValues(self, *args, **kwargs) -> 'CDSVector< float_type > const'
radii(self, *args, **kwargs) -> 'CDSVector< float_type >'
radiusScale(self, *args, **kwargs) -> 'float_type'
radiusType(self, *args, **kwargs) -> 'SolnScat::SolnScatRadiusType'
rho0(self, *args, **kwargs) -> 'float_type'
rhob(self, *args, **kwargs) -> 'float_type const'
rigidRegion(self, *args, **kwargs) -> 'SolnScatRigid const &'
selection(self, *args, **kwargs) -> 'AtomSel'
setAngleType(self, *args, **kwargs) -> 'void'
setBoundaryThickness(self, *args, **kwargs) -> 'void'
setCalcBoundary(self, *args, **kwargs) -> 'void'
setCalcType(self, *args, **kwargs) -> 'void'
setFormFactors(self, *args, **kwargs) -> 'void'
setFormScale(self, *args, **kwargs) -> 'void'
setFreeSurfVoxels(self, *args, **kwargs) -> 'void'
setGlobs(self, *args, **kwargs) -> 'void'
setIFormFactors(self, *args, **kwargs) -> 'void'
setInnerBoundaryFilename(self, *args, **kwargs) -> 'void'
setLMax(self, *args, **kwargs) -> 'void'
setNumAngles(self, *args, **kwargs) -> 'void'
setOuterBoundaryFilename(self, *args, **kwargs) -> 'void'
setQValues(self, *args, **kwargs) -> 'void'
setRadiusScale(self, *args, **kwargs) -> 'void'
setRadiusType(self, *args, **kwargs) -> 'void'
setRho0(self, *args, **kwargs) -> 'void'
setRhob(self, *args, **kwargs) -> 'void'
setSolventRadius(self, *args, **kwargs) -> 'void'
setUseGlobs(self, *args, **kwargs) -> 'void'
setVerbose(self, *args, **kwargs) -> 'void'
setVolumeScale(self, *args, **kwargs) -> 'void'
solventRadius(self, *args, **kwargs) -> 'float_type const'
surfVoxels(self, *args, **kwargs) -> 'CDSList< Voxel >'
useGlobs(self, *args, **kwargs) -> 'bool const'
verbose(self, *args, **kwargs) -> 'int const'
volumeScale(self, *args, **kwargs) -> 'float_type'

Static methods defined here:
__swig_destroy__ = delete_SolnScat(...)

Data descriptors defined here:
F

 
QValues

 
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
badVoxelFraction

 
boundaryF

 
boundaryI_

 
boundaryPoints

 
fGlob

 
fValues

 
posCenter

 
posGlob

 
rigidRegionModified

 
rigidRegions

 
surfVoxelsP

 
surfaceAreaFraction

 
thisown

 
The membership flag

Data and other attributes defined here:
CALC_CURRENT = 4
CALC_MULTIPOLE = 3
CALC_N = 1
CALC_N2 = 0
CALC_NONE = 6
CALC_REFERENCE = 5
CALC_UNIFORM = 2
RADIUS = 0
VOLUME = 1

 
class SolnScatPot_LetterClass(EnsemblePot)
    SolnScatPot_LetterClass(*args, **kwargs)
 

 
 
Method resolution order:
SolnScatPot_LetterClass
EnsemblePot
builtins.object

Methods defined here:
I_contrib(self, *args, **kwargs) -> 'CDSVector< float_type >'
__init__(self, *args, **kwargs)
Initialize self.  See help(type(self)) for accurate signature.
__repr__ = _swig_repr(self)
addRigidRegion(self, *args, **kwargs) -> 'void'
angleType(self, *args, **kwargs) -> 'XplorNIH::SpherePointType'
aveRadius(self, *args, **kwargs) -> 'float_type'
aveVolume(self, *args, **kwargs) -> 'float_type'
bg(self, *args, **kwargs) -> 'float_type'
boundaryI(self, *args, **kwargs) -> 'CDSVector< float_type >'
boundaryThickness(self, *args, **kwargs) -> 'float_type'
boundaryVol(self, *args, **kwargs) -> 'float_type'
calcBoundary(self, *args, **kwargs) -> 'bool'
calcBoundaryPoints(self, *args, **kwargs) -> 'void'
calcDerivs(self, *args, **kwargs) -> 'CDSVector< Vec3 >'
calcGlobCorrect_calcType(self, *args, **kwargs) -> 'void'
calcGlobCorrect_vec(self, *args, **kwargs) -> 'void'
calcType(self, *args, **kwargs) -> 'SolnScat::SolnScatCalcType'
calcd(self, *args, **kwargs) -> 'CDSVector< float_type >'
calcdScale(self, *args, **kwargs) -> 'float_type'
cmpType(self, *args, **kwargs) -> 'SolnScatPot::SolnScatCmpType'
convertFormFactors(self, *args, **kwargs) -> 'void'
energyMaybeDerivs0(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivs1(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivs2(self, *args, **kwargs) -> 'float_type'
expt(self, *args, **kwargs) -> 'CDSVector< float_type >'
exptScale(self, *args, **kwargs) -> 'float_type'
formFactors(self, *args, **kwargs) -> 'CDSVector< SolnScat::FormFactorMap > const &'
formScale(self, *args, **kwargs) -> 'float_type'
getF(self, *args, **kwargs) -> 'CDSMatrix< CDS::complex > const'
globCorrectCalcType(self, *args, **kwargs) -> 'SolnScat::SolnScatCalcType'
globCorrectMovement(self, *args, **kwargs) -> 'float_type const'
globCorrectNeedsUpdating(self, *args, **kwargs) -> 'bool const'
globCorrectReference(self, *args, **kwargs) -> 'rc_DerivedPot< SolnScatPot >'
globCorrectSteps(self, *args, **kwargs) -> 'int const'
iFormFactors(self, *args, **kwargs) -> 'CDSVector< SolnScat::FormFactorIMap > const &'
info(self, *args, **kwargs) -> 'String'
lMax(self, *args, **kwargs) -> 'int'
normalizeIndex(self, *args, **kwargs) -> 'int'
numAngles(self, *args, **kwargs) -> 'int'
numRestraints(self, *args, **kwargs) -> 'int'
pyXplorHelp(self, *args, **kwargs) -> 'String'
qValues(self, *args, **kwargs) -> 'CDSVector< float_type > const &'
qValuesExpt(self, *args, **kwargs) -> 'CDSVector< float_type > const &'
radii(self, *args, **kwargs) -> 'CDSVector< float_type >'
radiusScale(self, *args, **kwargs) -> 'float_type'
radiusType(self, *args, **kwargs) -> 'SolnScat::SolnScatRadiusType'
rho0(self, *args, **kwargs) -> 'float_type'
rhob(self, *args, **kwargs) -> 'float_type'
rms(self, *args, **kwargs) -> 'float_type'
selection(self, *args, **kwargs) -> 'AtomSel'
setAngleType(self, *args, **kwargs) -> 'void'
setBG(self, *args, **kwargs) -> 'void'
setBoundaryThickness(self, *args, **kwargs) -> 'void'
setCalcBoundary(self, *args, **kwargs) -> 'void'
setCalcType(self, *args, **kwargs) -> 'void'
setCmpType(self, *args, **kwargs) -> 'void'
setExpt(self, *args, **kwargs) -> 'void'
setFormFactors(self, *args, **kwargs) -> 'void'
setFormScale(self, *args, **kwargs) -> 'void'
setGlobCorrectCalcType(self, *args, **kwargs) -> 'void'
setGlobCorrectMovement(self, *args, **kwargs) -> 'void'
setGlobCorrectNeedsUpdating(self, *args, **kwargs) -> 'void'
setGlobCorrectReference(self, *args, **kwargs) -> 'void'
setGlobCorrectSteps(self, *args, **kwargs) -> 'void'
setGlobs(self, *args, **kwargs) -> 'void'
setIFormFactors(self, *args, **kwargs) -> 'void'
setLMax(self, *args, **kwargs) -> 'void'
setNormalizeIndex(self, *args, **kwargs) -> 'void'
setNumAngles(self, *args, **kwargs) -> 'void'
setQValues(self, *args, **kwargs) -> 'void'
setRadiusScale(self, *args, **kwargs) -> 'void'
setRadiusType(self, *args, **kwargs) -> 'void'
setRho0(self, *args, **kwargs) -> 'void'
setRhob(self, *args, **kwargs) -> 'void'
setSolventRadius(self, *args, **kwargs) -> 'void'
setUseGlobs(self, *args, **kwargs) -> 'void'
setVerbose(self, *args, **kwargs) -> 'void'
setVolumeScale(self, *args, **kwargs) -> 'void'
setWeights(self, *args, **kwargs) -> 'void'
simulation(self, *args) -> 'EnsembleSimulationBase const *'
solnScat(self, *args, **kwargs) -> 'SolnScat *'
solventRadius(self, *args, **kwargs) -> 'float_type'
splined(self, *args, **kwargs) -> 'CDSVector< float_type >'
updateDelta(self, *args, **kwargs) -> 'void'
updateGlobCorrect(self, *args, **kwargs) -> 'void'
useGlobs(self, *args, **kwargs) -> 'bool'
verbose(self, *args, **kwargs) -> 'int'
violations(self, *args, **kwargs) -> 'float_type'
volumeScale(self, *args, **kwargs) -> 'float_type'
weights(self, *args, **kwargs) -> 'CDSVector< float_type >'

Static methods defined here:
__swig_destroy__ = delete_SolnScatPot_LetterClass(...)

Data descriptors defined here:
dIs_dI

 
globCorrect

 
globCorrectCalcType_

 
thisown

 
The membership flag

Data and other attributes defined here:
CMP_LOG = 0
CMP_PLAIN = 1

Methods inherited from EnsemblePot:
addEnsWeights(self, *args, **kwargs) -> 'void'
calcEnergy(self, *args, **kwargs) -> 'float_type'
calcEnergyAndDerivs(self, *args, **kwargs) -> 'float_type'
calcWDerivs(self, *args, **kwargs) -> 'bool const'
clearEnsWeights(self, *args, **kwargs) -> 'void'
energyMaybeDerivs3(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivs4(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivsPost(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivsPre(self, *args, **kwargs) -> 'float_type'
ensWeight(self, *args, **kwargs) -> 'float_type'
ensWeights(self, *args, **kwargs) -> 'CDSList< float_type >'
ensWeightsInfo(self, *args, **kwargs) -> 'String'
getEnsWeights(self, *args, **kwargs) -> 'CDSList< VarEnsWeights > &'
setCalcWDerivs(self, *args, **kwargs) -> 'void'
setEnsWeights(self, *args, **kwargs) -> 'void'
setUseSimEnsWeights(self, *args, **kwargs) -> 'void'
updateEnsWeights(self, *args, **kwargs) -> 'void'
useSimEnsWeights(self, *args, **kwargs) -> 'bool const'

Data descriptors inherited from EnsemblePot:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)

 
class SolnScatRigid(builtins.object)
    SolnScatRigid(*args)
 

 
  Methods defined here:
__init__(self, *args)
Initialize self.  See help(type(self)) for accurate signature.
__repr__ = _swig_repr(self)
calcFContrib(self, *args, **kwargs) -> 'void'
calcFatQ(self, *args, **kwargs) -> 'complex'
derivContrib(self, *args, **kwargs) -> 'void'
update(self, *args, **kwargs) -> 'void'

Static methods defined here:
__swig_destroy__ = delete_SolnScatRigid(...)

Data descriptors defined here:
F

 
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
fitter

 
relAtomIndices

 
relGlobIndices

 
sphereInterp

 
thisown

 
The membership flag

 
class SolnScat_FracElec_SelIndex(builtins.object)
    SolnScat_FracElec_SelIndex(*args, **kwargs)
 

 
  Methods defined here:
__init__(self, *args, **kwargs)
Initialize self.  See help(type(self)) for accurate signature.
__repr__ = _swig_repr(self)

Static methods defined here:
__swig_destroy__ = delete_SolnScat_FracElec_SelIndex(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
fracElec

 
selIndex

 
thisown

 
The membership flag

 
class VarEnsWeights(builtins.object)
    VarEnsWeights(*args, **kwargs)
 

 
  Methods defined here:
__init__(self, *args, **kwargs)
Initialize self.  See help(type(self)) for accurate signature.
__repr__ = _swig_repr(self)

Static methods defined here:
__swig_destroy__ = delete_VarEnsWeights(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
ensWeights

 
mult

 
thisown

 
The membership flag

 
class Voxel(builtins.object)
    Voxel(*args, **kwargs)
 

 
  Methods defined here:
__init__(self, *args, **kwargs)
Initialize self.  See help(type(self)) for accurate signature.
__repr__ = _swig_repr(self)

Static methods defined here:
__swig_destroy__ = delete_Voxel(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag
v1

 
v1i

 
v2

 
v2i

 
v3

 
v3i

 
vol

 

 
class rc_EnsemblePot(builtins.object)
    rc_EnsemblePot(*args, **kwargs)
 

 
  Methods defined here:
__init__(self, *args, **kwargs)
Initialize self.  See help(type(self)) for accurate signature.
__repr__ = _swig_repr(self)

Static methods defined here:
__swig_destroy__ = delete_rc_EnsemblePot(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag

 
realSolnScatPot = class SolnScatPot(builtins.object)
    realSolnScatPot(*args)
 

 
  Methods defined here:
I_contrib(self, *args, **kwargs) -> 'CDSVector< float_type >'
__deref__(self, *args, **kwargs) -> 'SolnScatPot *'
__init__(self, *args)
__oldinit__ = __init__(self, *args, **kwargs)
__ref__(self, *args, **kwargs) -> 'SolnScatPot &'
__repr__ = _swig_repr(self)
addEnsWeights(self, *args, **kwargs) -> 'void'
addRigidRegion(self, *args, **kwargs) -> 'void'
angleType(self, *args, **kwargs) -> 'XplorNIH::SpherePointType'
aveRadius(self, *args, **kwargs) -> 'float_type'
aveVolume(self, *args, **kwargs) -> 'float_type'
bg(self, *args, **kwargs) -> 'float_type'
boundaryI(self, *args, **kwargs) -> 'CDSVector< float_type >'
boundaryThickness(self, *args, **kwargs) -> 'float_type'
boundaryVol(self, *args, **kwargs) -> 'float_type'
calcBoundary(self, *args, **kwargs) -> 'bool'
calcBoundaryPoints(self, *args, **kwargs) -> 'void'
calcDerivs(self, *args, **kwargs) -> 'CDSVector< Vec3 >'
calcEnergy(self, *args, **kwargs) -> 'float_type'
calcEnergyAndDerivs(self, *args, **kwargs) -> 'float_type'
calcGlobCorrect(s, arg='current')
swig can't handle the enumtypemap with this overridden method:
we must do it manually
calcGlobCorrect_calcType(self, *args, **kwargs) -> 'void'
calcGlobCorrect_vec(self, *args, **kwargs) -> 'void'
calcType(self, *args, **kwargs) -> 'SolnScat::SolnScatCalcType'
calcWDerivs(self, *args, **kwargs) -> 'bool const'
calcd(self, *args, **kwargs) -> 'CDSVector< float_type >'
calcdScale(self, *args, **kwargs) -> 'float_type'
clearEnsWeights(self, *args, **kwargs) -> 'void'
cmpType(self, *args, **kwargs) -> 'SolnScatPot::SolnScatCmpType'
convertFormFactors(self, *args, **kwargs) -> 'void'
decrRefCnt(self, *args, **kwargs) -> 'void'
energyMaybeDerivs0(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivs1(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivs2(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivs3(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivs4(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivsPost(self, *args, **kwargs) -> 'float_type'
energyMaybeDerivsPre(self, *args, **kwargs) -> 'float_type'
ensWeight(self, *args, **kwargs) -> 'float_type'
ensWeights(self, *args, **kwargs) -> 'CDSList< float_type >'
ensWeightsInfo(self, *args, **kwargs) -> 'String'
expt(self, *args, **kwargs) -> 'CDSVector< float_type >'
exptScale(self, *args, **kwargs) -> 'float_type'
formFactors(self, *args, **kwargs) -> 'CDSVector< SolnScat::FormFactorMap > const &'
formScale(self, *args, **kwargs) -> 'float_type'
getEnsWeights(self, *args, **kwargs) -> 'CDSList< VarEnsWeights > &'
getF(self, *args, **kwargs) -> 'CDSMatrix< CDS::complex > const'
get_weights(self, *args, **kwargs) -> 'CDSVector< float_type >'
globCorrectCalcType(self, *args, **kwargs) -> 'SolnScat::SolnScatCalcType'
globCorrectMovement(self, *args, **kwargs) -> 'float_type const'
globCorrectNeedsUpdating(self, *args, **kwargs) -> 'bool const'
globCorrectReference(self, *args, **kwargs) -> 'rc_DerivedPot< SolnScatPot >'
globCorrectSteps(self, *args, **kwargs) -> 'int const'
iFormFactors(self, *args, **kwargs) -> 'CDSVector< SolnScat::FormFactorIMap > const &'
incrRefCnt(self, *args, **kwargs) -> 'void'
info(self, *args, **kwargs) -> 'String'
instanceData(self, *args, **kwargs) -> 'PyObject *'
instanceName(self, *args, **kwargs) -> 'char const *'
lMax(self, *args, **kwargs) -> 'int'
normalizeIndex(self, *args, **kwargs) -> 'int'
numAngles(self, *args, **kwargs) -> 'int'
numRestraints(self, *args, **kwargs) -> 'int'
potName(self, *args, **kwargs) -> 'char const *'
pyXplorHelp(self, *args, **kwargs) -> 'String'
qValues(self, *args, **kwargs) -> 'CDSVector< float_type > const &'
qValuesExpt(self, *args, **kwargs) -> 'CDSVector< float_type > const &'
radii(self, *args, **kwargs) -> 'CDSVector< float_type >'
radiusScale(self, *args, **kwargs) -> 'float_type'
radiusType(self, *args, **kwargs) -> 'SolnScat::SolnScatRadiusType'
refCnt(self, *args, **kwargs) -> 'int'
registerInstanceData(self, *args, **kwargs) -> 'void'
registerTo(self, *args, **kwargs) -> 'void'
resetInstanceName(self, *args, **kwargs) -> 'void'
resetPotName(self, *args, **kwargs) -> 'void'
rho0(self, *args, **kwargs) -> 'float_type'
rhob(self, *args, **kwargs) -> 'float_type'
rms(self, *args, **kwargs) -> 'float_type'
scale(self, *args, **kwargs) -> 'float_type const'
selection(self, *args, **kwargs) -> 'AtomSel'
setAngleType(self, *args, **kwargs) -> 'void'
setBG(self, *args, **kwargs) -> 'void'
setBoundaryThickness(self, *args, **kwargs) -> 'void'
setCalcBoundary(self, *args, **kwargs) -> 'void'
setCalcType(self, *args, **kwargs) -> 'void'
setCalcWDerivs(self, *args, **kwargs) -> 'void'
setCmpType(self, *args, **kwargs) -> 'void'
setEnsWeights(self, *args, **kwargs) -> 'void'
setExpt(self, *args, **kwargs) -> 'void'
setFormFactors(self, *args, **kwargs) -> 'void'
setFormScale(self, *args, **kwargs) -> 'void'
setGlobCorrectCalcType(self, *args, **kwargs) -> 'void'
setGlobCorrectMovement(self, *args, **kwargs) -> 'void'
setGlobCorrectNeedsUpdating(self, *args, **kwargs) -> 'void'
setGlobCorrectReference(self, *args, **kwargs) -> 'void'
setGlobCorrectSteps(self, *args, **kwargs) -> 'void'
setGlobs(self, *args, **kwargs) -> 'void'
setIFormFactors(self, *args, **kwargs) -> 'void'
setLMax(self, *args, **kwargs) -> 'void'
setNormalizeIndex(self, *args, **kwargs) -> 'void'
setNumAngles(self, *args, **kwargs) -> 'void'
setQValues(self, *args, **kwargs) -> 'void'
setRadiusScale(self, *args, **kwargs) -> 'void'
setRadiusType(self, *args, **kwargs) -> 'void'
setRho0(self, *args, **kwargs) -> 'void'
setRhob(self, *args, **kwargs) -> 'void'
setScale(self, *args, **kwargs) -> 'void'
setSolventRadius(self, *args, **kwargs) -> 'void'
setThreshold(self, *args, **kwargs) -> 'void'
setUseGlobs(self, *args, **kwargs) -> 'void'
setUseSimEnsWeights(self, *args, **kwargs) -> 'void'
setVerbose(self, *args, **kwargs) -> 'void'
setVolumeScale(self, *args, **kwargs) -> 'void'
setWeights(self, *args, **kwargs) -> 'void'
set_weights(self, *args, **kwargs) -> 'void'
simulation(self, *args) -> 'EnsembleSimulationBase const *'
solnScat(self, *args, **kwargs) -> 'SolnScat *'
solventRadius(self, *args, **kwargs) -> 'float_type'
splined(self, *args, **kwargs) -> 'CDSVector< float_type >'
threshold(self, *args, **kwargs) -> 'float_type const'
unRegister(self, *args, **kwargs) -> 'void'
updateDelta(self, *args, **kwargs) -> 'void'
updateEnsWeights(self, *args, **kwargs) -> 'void'
updateGlobCorrect(self, *args, **kwargs) -> 'void'
updateValues(self, *args, **kwargs) -> 'void'
useGlobs(self, *args, **kwargs) -> 'bool'
useSimEnsWeights(self, *args, **kwargs) -> 'bool const'
verbose(self, *args, **kwargs) -> 'int'
violations(self, *args, **kwargs) -> 'float_type'
volumeScale(self, *args, **kwargs) -> 'float_type'
weights(self, *args, **kwargs) -> 'CDSVector< float_type >'

Static methods defined here:
__swig_destroy__ = delete_SolnScatPot(...)

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
dIs_dI

 
globCorrect

 
globCorrectCalcType_

 
instanceDataCleanup

 
instanceDataCreate

 
instanceData_

 
modified

 
registeredSimulations

 
thisown

 
The membership flag

 
Functions
       
SolnScatPot(*args)
calcIfromF(*args, **kwargs) -> 'CDSVector< float_type,0,CDS::DefaultAlloc,long,int >'
calc_dIdrhob(*args, **kwargs) -> 'CDSVector< float_type,0,CDS::DefaultAlloc,long,int >'
eigen(*args) -> 'CDSList< EigenPair< double >,10 >'
fast_max(*args, **kwargs) -> 'CDSMatrix< CDS::complex >::ElementType'
fast_min(*args, **kwargs) -> 'CDSMatrix< CDS::complex >::ElementType'
fromPy(*args, **kwargs) -> 'void'
getCol(*args, **kwargs) -> 'CDSVector< CDSMatrix< CDS::complex >::ElementType >'
getRow(*args, **kwargs) -> 'CDSVector< CDSMatrix< CDS::complex >::ElementType >'
pairDistribution(*args, **kwargs) -> 'CDSList< FloatPair,10 >'
pointsOnSphere(*args, **kwargs) -> 'CDSVector< Vec3 >'
pyXplorHelp(*args) -> 'String'
transpose(*args, **kwargs) -> 'CDSMatrix< CDS::complex >'
vec_outerProd(*args, **kwargs) -> 'CDSMatrix< CDS::complex >'
volContrib(*args, **kwargs) -> 'CDSMatrix< CDS::complex >'

 
Data
        ANGLE_FIBONACCI = 1
ANGLE_RANDOM = 0
ANGLE_SPIRAL = 2
__CDSVector_hh__ = 1