[table of contents] [master index] [comments] [modules] [programs] [variables] [types] [procedures]
MODULE
sweep_mod
PURPOSE
Solve SCF repeatedly for a range of parameters, using
zeroth, first, or second order continuation of solutions
COMMENT
Module uses 1st or 2nd order contination of scf equations
at equally spaced values of a continuation parameter s along
a line in parameter space along which each of the chemistry
input parameters kuhn, block_length, chi or Temperature
(depending upon the value of chi_fmt), and phi or mu
(depending upon the ensemble) has a constant increment per
unit change in s. For solutions of a specified sequence of
unit cell parameters, with domain = F, (rather than iteration
of each set of chemical parameters to stress free solution,
with domain = T), the elements of cell_param also have a
constant increment per unit change in s. A single parameter
may be varied by setting all but one of the increments to zero.
The required increments are read from file by subroutine
input_increments and are stored in private module variables
d_kuhn, d_block_length, d_chi or d_Temperature, etc.
The order of continuation is set at compile time by the value
of the private module parameter N_hist, which gives the maximum
number of solutions stored in a history stack. Compile with
N_hist = 3 (current solution and two previous solutions) for
second order continuation, and N_hist = 2 (current and one
previous solution) for first order continuation.
SOURCE
module sweep_mod use const_mod use unit_cell_mod, only : N_cell_param use chemistry_mod implicit none private ! Public procedures public :: input_increments ! Read nonzero increments from file public :: output_increments ! Write increments to file public :: increment_parameters ! Increment parameters public :: history_setup ! Allocate history stack arrays public :: update_history ! Push new solution onto history stack public :: continuation ! Predict next solution
SUBROUTINE
input_increments(i_unit,fmt_flag,domain)
PURPOSE
Allocate arrays for chemistry increment variables
Read values of increment module variables from unit i_unit
NOTE
Default imput unit in module io_mod is set to i_unit
SOURCE
subroutine input_increments(i_unit,fmt_flag,domain) use io_mod integer, intent(IN) :: i_unit ! input unit # character(1), intent(IN) :: fmt_flag ! F = formatted, U = unformatted logical, intent(IN) :: domain ! true if domain iteration
SUBROUTINE
output_increments(o_unit,fmt_flag,domain)
PURPOSE
Write values of increments to unit o_unit
in format required by input_increments
NOTE
Default output unit in module io_mod set to o_unit
SOURCE
subroutine output_increments(o_unit,fmt_flag,domain) use io_mod integer, intent(IN) :: o_unit ! input unit # character(1), intent(IN) :: fmt_flag ! F = formatted, U = unformatted logical, intent(IN) :: domain ! true if domain iteration
SUBROUTINE
increment_parameters(step,domain,cell_param)
PURPOSE
Increment chemisty parameters
if (.not.domain), also increment cell_param
ARGUMENTS
step - change in continuation parameter s
domain - F -> specified unit cell, T -> variable unit cell
cell_param - cell parameters
SOURCE
subroutine increment_parameters(step,domain,cell_param) real(long) :: step logical :: domain real(long) :: cell_param(:)
SUBROUTINE
history_setup
PURPOSE
Allocates private history stack module variables:
s_hist(N_hist) = prior values of s
omega_hist(N_monomer,N_star,N_hist) = prior omega fields
cell_param_hist(6, N_hist) = prior cell parameters
SOURCE
subroutine history_setup use basis_mod, only : N_star
SUBROUTINE
update_history(s,omega,cell_param,domain)
PURPOSE
Push new values of s, omega, & cell_param onto history stacks.
On output:
s_hist(1) = s
omega_hist(:,:,1) = omega
cell_param_hist(:,1) = cell_param
The cell_param stack is updated only if domain = T
SOURCE
subroutine update_history(s,omega,cell_param,domain) use basis_mod, only : N_star real(long), intent(IN) :: s real(long), intent(IN) :: omega(:,:) real(long), intent(IN) :: cell_param(:) logical, intent(IN) :: domain
SUBROUTINE
continuation(step,domain,omega,cell_param)
PURPOSE
Predict new values of omega and (if domain = T) cell_param
for a change in s of size step. Uses the history stack of
previous values of omega, cell parameters, and s.
Order of continuation depends on the value of the private
variable depth_hist = current depth of stack:
If depth_hist = 1 -> zeroth order continuation
If depth_hist = 2 -> first order continuation
If depth_hist > 2 -> second order continuation
SOURCE
subroutine continuation(step,domain,omega,cell_param) real(long), intent(IN) :: step ! change in s logical, intent(IN) :: domain ! adjust cell_param real(long), intent(INOUT) :: omega(:,:) ! chemical potential field real(long), intent(INOUT) :: cell_param(:)! unit cell parameters