Polymer Field Theory: Liquid Crystalline Solutions

We define a system Hamiltonian [Spakowitz2003] of a collection of polymer chains and solvent molecules

(1)\[\beta \mathcal{H} = \sum_{i=1}^{n_{p}} \frac{\epsilon}{2} \int_{0}^{L} \left( \frac{\partial \vec{u}_{i}}{\partial s} \right)^{\! \! 2} d s + \chi \int d \vec{r} \, \hat{\phi}_{s}(\vec{r} \,) \hat{\phi}_{p}(\vec{r} \,) - \frac{a}{2} \int d \vec{r} \, \hat{\mathbf{S}}(\vec{r} \,):\hat{\mathbf{S}}(\vec{r} \,)\]

\(\hat{\phi}_{s}\) and \(\hat{\phi}_{p}\) are the local volume fractions of the solvent and polymer molecules, and \(\hat{\mathbf{S}}\) is the tensorial nematic order parameter

\[\begin{split}\hat{\phi}_{s}(\vec{r} \,) & = & v \sum_{j=1}^{n_{s}} \delta (\vec{r}-\vec{r}_{j}) \\ \hat{\phi}_{p}(\vec{r} \,) & = & A \sum_{i=1}^{n_{p}} \int_{0}^{L} d s \delta (\vec{r}-\vec{r}_{i}(s)) \\ \hat{\mathbf{S}} (\vec{r} \,) & = & A \sum_{i=1}^{n_{p}} \int_{0}^{L} d s \delta (\vec{r}-\vec{r}_{i}(s)) %\nonumber \\ \left( \vec{u}_{i}(s) \vec{u}_{i}(s) - \frac{1}{3} \mathbf{I} \right)\end{split}\]

The Flory-Huggins parameter \(\chi\) gives the strength of polymer-solvent interactions. The Maier-Saupe parameter \(a\) gives the strength of aligning interactions that lead to liquid crystal phase phenomena.

The thermodynamic behavior is found by finding the grand canonical partition function

(2)\[\begin{split}\Xi & = & \exp [-\beta G(T,V,\mu)] \nonumber \\ & = &\sum_{n_{s},n_{p}=0}^{\infty} \frac{1}{n_{s}! n_{p}!} \frac{e^{\beta \mu L A n_{p}}}{v^{n_{s}} (L A)^{n_{p}}} \int \prod_{j=1}^{n_{s}} d \vec{r}_{j} \int \prod_{i=1}^{n_{p}} \mathcal{D}[\vec{r}_{i}] \nonumber \\ & & \hspace{0.1in} \times \prod_{\vec{r}} \delta (\hat{\phi}_{s}+\hat{\phi}_{p}-1) \prod_{s}\delta (\vec{u}_{i} - \partial_{s} \vec{r}_{i}) e^{-\beta \mathcal{H}[\hat{\phi}_{s},\hat{\phi}_{p},\hat{\mathbf{S}}]}\end{split}\]

where \(\prod_{\vec{r}} \delta (\hat{\phi}_{s}+\hat{\phi}_{p}-1)\) accounts for the incompressibility constraint \(\hat{\phi}_{s}+\hat{\phi}_{p} = 1\) at all locations, and \(\prod_{s}\delta (\vec{u}_{i} - \partial_{s} \vec{r}_{i})\) denotes the fixed chain length constraint \(\vec{u}_{i}=\partial_{s} \vec{r}_{i}\). This cannot be solved exactly due to the many-body interactions that are implicit within the Hamiltonian.

To make progress, we use field-theoretical techniques~cite{kn:fredreview,kn:matsen} to transform the many-chain problem into a single-chain problem in fluctuating effective potential fields. The field-theoretic representation emerges from a series of identity transformations that leave the theory unaffected.

First, we introduce collective variables \(\phi_{s}\), \(\phi_{p}\), and \(\mathbf{S}\) into the grand canonical partition function, giving

\[\begin{split}\Xi & = & \int \! \! \mathcal{D} \phi_{s} \int \! \! \mathcal{D} \phi_{p} \int \! \! \mathcal{D} \mathbf{S} \! \! \! \sum_{n_{s},n_{p}=0}^{\infty} \frac{1}{n_{s}! n_{p}!} \frac{e^{\beta \mu L A n_{p}}}{v^{n_{s}} (L A)^{n_{p}}} \int \prod_{j=1}^{n_{s}} d \vec{r}_{j} \nonumber \\ & & \hspace{0.1in} \times \int \prod_{i=1}^{n_{p}} \mathcal{D}[\vec{r}_{i}] \prod_{\vec{r}} \delta (\phi_{s} - \hat{\phi}_{s}) \prod_{\vec{r}} \delta (\phi_{p} - \hat{\phi}_{p}) \prod_{\vec{r}} \delta (\mathbf{S} - \hat{\mathbf{S}}) \nonumber \\ & & \hspace{0.1in} \times \prod_{\vec{r}} \delta (\phi_{s}+\phi_{p}-1) \prod_{s}\delta (\vec{u}_{i} - \partial_{s} \vec{r}_{i}) e^{-\beta \mathcal{H}[\phi_{s},\phi_{p},\mathbf{S}]}\end{split}\]

Second, rewrite the delta functions using Fourier representations

\[\begin{split}\prod_{\vec{r}} \delta (\phi_{s} - \hat{\phi}_{s}) = \int \! \! \mathcal{D} W_{s} \exp \left[ i \int \! \! d \vec{r} W_{s}(\vec{r}) \left( \phi_{s}(\vec{r}) - \hat{\phi}_{s}(\vec{r}) \right) \right] \\ \prod_{\vec{r}} \delta (\phi_{p} - \hat{\phi}_{p}) = \int \! \! \mathcal{D} W_{p} \exp \left[ i \int \! \! d \vec{r} W_{p}(\vec{r}) \left( \phi_{p}(\vec{r}) - \hat{\phi}_{p}(\vec{r}) \right) \right] \\ \prod_{\vec{r}} \delta (\mathbf{S} - \hat{\mathbf{S}}) = \int \! \! \mathcal{D} \mathbf{\lambda} \exp \left[ i \int \! \! d \vec{r} \mathbf{\lambda}(\vec{r}) : \left( \mathbf{S}(\vec{r}) - \hat{\mathbf{S}}(\vec{r}) \right) \right]\end{split}\]

The grand canonical partition function is now written as

(3)\[\begin{split}\Xi & = & \int \mathcal{D} \phi_{p} \mathcal{D} W_{s} \mathcal{D} W_{p} \mathcal{D} \mathbf{S} \mathcal{D} \mathbf{\lambda} \exp \left\{ i \int d \vec{r} \, W_{s} (1-\phi_{p}) + \right. \nonumber \\ & & i \int d \vec{r} \, \left[ W_{p} \phi_{p} + \mathbf{\lambda}:\mathbf{S} \right] - \chi \int d \vec{r} \, \phi_{p} (1-\phi_{p}) + \nonumber \\ & & \left. \frac{a}{2} \int d \vec{r} \, \mathbf{S}:\mathbf{S} + \frac{1}{v} z_{s}(W_{s}) + \frac{e^{\beta \mu L A}}{L A} z_{p}(W_{p}, \mathbf{\lambda}) \right\}\end{split}\]

Here, we have defined the single-solvent-molecule partition function and single-polymer-chain partition function, given by

(4)\[z_{s}(W_{s}) = \int d \vec{r} \, e^{-i v W_{s}(\vec{r} \,)}\]

and

(5)\[\begin{split}& & z_{p}(W_{p},\mathbf{\lambda}) = \int \mathcal{D}[\vec{r} \,(s)] \exp \left\{- \frac{\epsilon}{2} \int_{0}^{L} \left( \frac{\partial \vec{u}}{\partial s} \right)^{\! \! 2} d s - \right. \nonumber \\ & & i A \int_{0}^{L} d s \left[ W_{p}(\vec{r} \,(s)) + \left. \mathbf{\lambda}(\vec{r} \,(s)) \! : \! \! \left(\vec{u} \vec{u} - \frac{1}{3} \mathbf{I}\right) \right] \right\}\end{split}\]

This form still cannot be solved exactly, but we can find a saddle-point approximation that essentially gives the maximum term in the functional integrals. To do this, we write the free energy functional

\[\begin{split}& & - \beta G = i \int d \vec{r} \, W_{s} (1-\phi_{p}) + i \int d \vec{r} \, \left[ W_{p} \phi_{p} + \mathbf{\lambda}:\mathbf{S} \right] \\ & & - \chi \int d \vec{r} \, \phi_{p} (1-\phi_{p}) + \frac{a}{2} \int d \vec{r} \, \mathbf{S}:\mathbf{S} + \frac{1}{v} z_{s}(W_{s}) + \frac{e^{\beta \mu L A}}{L A} z_{p}(W_{p}, \mathbf{\lambda}) \nonumber\end{split}\]

which will be evaluated at values of \(\phi_{s}\), \(\phi_{p}\), \(\mathbf{S}\), \(W_{s}\), \(W_{p}\), and \(\mathbf{\lambda}\) that maximize \(-\beta G\) (or free energy minimum).

We find the first variation of the free energy to be

\[\begin{split}& & - \delta \beta G = - i \int \! \! d \vec{r} W_{s} \delta \phi_{p} +i \int \! \! d \vec{r} (1- \phi_{p}) \delta W_{s} + i \int \! \! d \vec{r} W_{p} \delta \phi_{p} + i \int \! \! d \vec{r} \phi_{p} \delta W_{p} \nonumber \\ & & + i \int \! \! d \vec{r} \mathbf{\lambda} : \delta \mathbf{S} + i \int \! \! d \vec{r} \mathbf{S} : \delta \mathbf{\lambda} - \chi \int \! \! d \vec{r} (1 - 2 \phi_{p}) \delta \phi_{p} + a \int \! \! d \vec{r} \mathbf{S} : \delta \mathbf{S} \nonumber \\ & & +\frac{1}{v} \delta z_{s}(W_{s}) + \frac{e^{\beta \mu L A}}{L A} \delta z_{p}(W_{p}, \mathbf{\lambda})\end{split}\]

where

\[\begin{split}& & \delta z_{s}(W_{s}) = - i v \int \! \! d \vec{r} \exp \left[ - i v W_{s}(\vec{r}) \right] \delta W_{s} \\ & & \delta z_{p}(W_{p}, \mathbf{\lambda}) = \int \! \! d \vec{r} \int \mathcal{D}[\vec{r} \,(s)] \int_{0}^{L} \! \! ds \left[ -iA \delta (\vec{r} - \vec{r}(s)) \delta W_{p}(\vec{r}) - i A \delta (\vec{r} - \vec{r}(s)) \left(\vec{u} \vec{u} - \frac{1}{3} \mathbf{I}\right) : \delta \mathbf{\lambda} \right] \nonumber \\ & & \exp \left\{- \frac{\epsilon}{2} \int_{0}^{L} \left( \frac{\partial \vec{u}}{\partial s} \right)^{\! \! 2} d s - i A \int_{0}^{L} d s \left[ W_{p}(\vec{r} \,(s)) + \mathbf{\lambda}(\vec{r} \,(s)) \! : \! \! \left(\vec{u} \vec{u} - \frac{1}{3} \mathbf{I}\right) \right] \right\}\end{split}\]

Setting \(\delta \beta G\) to zero gives the self-consistent field equations:

\[\begin{split}1-\phi_{p} & = & e^{-v w_{s}}, \\ \phi_{p} & = & - \frac{e^{\beta \mu L A}}{L A} \frac{\delta z_{p}}{\delta w_{p}}, \\ w_{p} - w_{s} & = & \chi (1 - 2 \phi_{p}), \\ \mathbf{h} & = & - a \mathbf{S}, \label{saddles} \\ \mathbf{S} & = & - \frac{ e^{\beta \mu L A}}{L A} \frac{\delta z_{p}}{\delta \mathbf{h}},\end{split}\]

where \(w_{s} = i W_{s}\), \(w_{p} = i W_{p}\), and \(\mathbf{h}=i \mathbf{\lambda}\).

Now consider the special case of having a homogeneous solution \(\phi_{s}(\vec{r}) = \phi_{s}\), \(\phi_{p}(\vec{r}) = \phi_{p}\), and nematic order parameter \(\mathbf{S}(\vec{r}) = S_{0} \left( \hat{z}\hat{z} - \frac{1}{3} \mathbf{I} \right)\). Define a normalized scalar order parameter \(m\) by

(6)\[m = \frac{S_{0}}{\phi_{p}} = \frac{1}{L} \int_{0}^{L} d s \left< \left( \frac{3}{2} u_{z}^{2}(s) - \frac{1}{2} \right) \right>,\]

where \(\left< \ldots \right>\) indicates an average with respect to the single-chain, self-consistent-field Hamiltonian,

(7)\[\beta \mathcal{H}_{0} = \frac{\epsilon}{2} \int_{0}^{L} \left( \frac{\partial \vec{u}}{\partial s} \right)^{2} d s - a \phi_{p} m A \int_{0}^{L} d s \left( u_{z}^{2}(s) - \frac{1}{3} \right).\]

The grand potential density \(g = G/V\) is equivalently expressed as the osmotic pressure defined through \(p=-[g(\phi_p)-g(0)]\) and given by

(8)\[\begin{split}\beta p & = & -\frac{1}{v} \log (1-\phi_{p}) - \chi \phi_{p}^{2} -\frac{1}{3} a \phi_{p}^{2} m^{2} \nonumber \\ & & \hspace{0.5in} - \frac{\phi_{p}}{v} + \frac{\phi_{p}}{L A}\end{split}\]

The chemical potential \(\mu\) is given by

(9)\[\begin{split}\beta \mu & = & \frac{1}{L A} \log \phi_{p} - \frac{1}{v} \log (1-\phi_{p})+ \nonumber \\ & & \hspace{0.5in} \chi (1-2 \phi_{p}) - \frac{1}{L A} \log q\end{split}\]

where the single-polymer molecule orientation partition function \(q\) is calculated from the self-consistent field Hamiltonian, Eq. (7). The Helmholtz free energy density (up to an additive constant) for a single, homogeneous phase is

(10)\[\begin{split}\beta f & = & -\beta p + \beta \mu \phi_{p} \nonumber \\ & = & \frac{\phi_{p}}{L A} \log \phi_{p} + \frac{(1-\phi_{p})}{v} \log (1-\phi_{p}) +\chi \phi_{p} (1-\phi_{p}) \nonumber \\ & & + \frac{a}{3} \phi_{p}^{2} m^{2} +\frac{\phi_{p}}{v} - \frac{\phi_{p}}{L A} - \frac{\phi_{p}}{L A} \log q\end{split}\]

This agrees with the Flory-Huggins theory when \(a=0\) and \(q\) is a constant (with unimportant factors that are linear in \(\phi_{p}\)).

Functions contained with the ‘wlc_lcpoly’ module

wlc_lcpoly

Module containing functions for evaluating the self-consistent field theory for polymer solutions, including nematic liquid crystallinity

wlcstat.wlc_lcpoly.calc_int_mag(length_kuhn, poles_vec)[source]

Evaluate the magnitude of the integral for a list of poles (including repeats). This algorithm includes cases for single, double, and triple poles (as needed in evaluation of correlation functions)

Parameters
  • length_kuhn (float) – The length of the chain in Kuhn lengths

  • poles_vec (float (array)) – Array of poles

Returns

int_mag – Value of the integral over chain length for the five poles

Return type

float

wlcstat.wlc_lcpoly.elastic_lcpoly(length_kuhn, lam, alpha_max=25, l_cont_frac_max=200)[source]

Calculate the Frank elastic constants for a polymer liquid crystal solution

Parameters
  • length_kuhn (float) – The length of the chain in Kuhn lengths

  • lam (float) – The value of the quadrupole field \(\lambda\)

  • alpha_max (int) – Maximum number of poles evaluated (default 50)

  • l_cont_frac_max (int) – Number of levels included in the evaluation of residues from continued fraction

Returns

  • q_val (float) – Single-chain partition function

  • m_val (float) – Nematic order parameter

  • y21_y21 (float) – y21-y21 correlation function

  • y21_ux_ux_y21 (float) – y21-y21 correlation function with x-x end-to-end distance squared

  • y21_uy_uy_y21 (float) – y21-y21 correlation function with y-y end-to-end distance squared

  • y21_uz_uz_y21 (float) – y21-y21 correlation function with z-z end-to-end distance squared

wlcstat.wlc_lcpoly.elastic_rr(length_kuhn, lam)[source]

Calculate the Frank elastic constants for a polymer liquid crystal solution

Parameters
  • length_kuhn (float) – The length of the chain in Kuhn lengths

  • lam (float) – The value of the quadrupole field \(\lambda\)

Returns

  • q_val (float) – Single-chain partition function

  • m_val (float) – Nematic order parameter

  • y21_y21 (float) – y21-y21 correlation function

  • y21_ux_ux_y21 (float) – y21-y21 correlation function with x-x end-to-end distance squared

  • y21_uy_uy_y21 (float) – y21-y21 correlation function with y-y end-to-end distance squared

  • y21_uz_uz_y21 (float) – y21-y21 correlation function with z-z end-to-end distance squared

wlcstat.wlc_lcpoly.eval_a_l_m(l, m)[source]

eval_a_l_m - Evaluate the coefficient from a ladder operation \(cos \theta Y_{\lambda;\mu}\) on the spherical harmonic

Parameters
  • l (int (array)) – The angular kinetic energy quantum index of the spherical harmonic \(Y_{\lambda;\mu}\)

  • m (int) – The angular kinetic energy quantum index of the spherical harmonic \(Y_{\lambda;\mu}\)

Notes

See Mehraeen, et al, Phys. Rev. E, 77, 061803 (2008). (Ref [Mehraeen2008]) and Arfken (1999) (Ref [Arfken1999])

wlcstat.wlc_lcpoly.eval_b_l_m(l, m)[source]

eval_a_l_m - Evaluate the coefficient from a ladder operation \(e^{i \phi} sin \theta Y_{\lambda;\mu}\) on the spherical harmonic

Parameters
  • l (int (array)) – The angular kinetic energy quantum index of the spherical harmonic \(Y_{\lambda;\mu}\)

  • m (int) – The angular kinetic energy quantum index of the spherical harmonic \(Y_{\lambda;\mu}\)

Notes

See Mehraeen, et al, Phys. Rev. E, 77, 061803 (2008). (Ref [Mehraeen2008]) and Arfken (1999) (Ref [Arfken1999])

wlcstat.wlc_lcpoly.eval_poles_lcpoly(lam, m=0, alpha_max=25)[source]

eval_poles_lcpoly - Evaluate the poles for given \(\lambda\) and \(\mu\) using the matrix method for intermediate \(K\)

Parameters
  • lam (float) – The value of the nematic field \(\lambda\)

  • m (int) – Value of the mu parameter

  • alpha_max (int) – Maximum number of poles evaluated (default 25)

Returns

poles – Evaluated poles for the given \(\lambda\) and \(\mu\)

Return type

float

Notes

wlcstat.wlc_lcpoly.eval_residues_lcpoly(lam, m, poles, l_zero_only=True, l_max=25, alpha_max=25, l_cont_frac_max=50)[source]

eval_residues_lcpoly - Evaluate the residues for the Green’s function of a nematic polymer

Parameters
  • lam (float) – The value of the nematic field \(\lambda\)

  • m (int) – Value of the mu parameter

  • poles (complex float) – Evaluated poles for the given \(K\) and \(\mu\)

  • l_zero_only (boolean) – Indicates whether the residues will be evaluated over the range of \(\lambda\) and \(lambda_{0}\)

  • l_max (int) – Maximum lambda value evaluated

  • alpha_max (int) – Maximum number of poles evaluated (default 25)

  • l_cont_frac_max (int) – Maximum \(\lambda\) value in the continued fraction evaluation

Returns

residues – Evaluated residues for the given \(K\) and \(\mu\)

Return type

complex float

Notes

See [Mehraeen2008] for intermediate-k algorithms

wlcstat.wlc_lcpoly.m_lcpoly(length_kuhn, lam, alpha_max=25, l_cont_frac_max=50)[source]

Calculate the single-polymer partition function

Parameters
  • length_kuhn (float) – The length of the chain in Kuhn lengths

  • lam (float) – The value of the quadrupole field \(\lambda\)

  • alpha_max (int) – Maximum number of poles evaluated (default 50)

Returns

m_val – Evaluated order parameter

Return type

float

wlcstat.wlc_lcpoly.q_lcpoly(length_kuhn, lam, alpha_max=25)[source]

Calculate the single-polymer partition function

Parameters
  • length_kuhn (float) – The length of the chain in Kuhn lengths

  • lam (float) – The value of the quadrupole field \(\lambda\)

  • alpha_max (int) – Maximum number of poles evaluated (default 50)

Returns

q_val – Evaluated single-polymer partition function

Return type

float

wlcstat.wlc_lcpoly.r_2_lcpoly(length_kuhn, lam, alpha_max=25, l_cont_frac_max=50)[source]

Calculate the mean-square end-to-end distance for a liquid crystal polymer

Parameters
  • length_kuhn (float) – The length of the chain in Kuhn lengths

  • lam (float) – The value of the quadrupole field \(\lambda\)

  • alpha_max (int) – Maximum number of poles evaluated (default 50)

Returns

  • r_2_par (float) – mean-square distance in the parallel direction

  • r_2_perp (float) – mean-square distance in the perpendicular direction

  • xi_par (float) – parallel correlation length

  • xi_perp (float) – perpendicular correlation length

Example usage of ‘q_lcpoly’

We reproduce Fig. 6 from [Spakowitz2003] to demonstrate the use of ‘q_lcpoly’. We show free energy relative to the isotropic state \(\Delta f\) versus the order parameter \(m\) for a thermotropic liquid-crystalline polymer system with \(N = L/(2lp)=3.3333\) and \(2 l_{p} \kappa = 20.7582\) (blue), \(2 l_{p} \kappa = 21.0606\) (orange), and \(2 l_{p} \kappa = 23.6844\) (green).

(Source code)