NLSL
ordrpr.f90
Go to the documentation of this file.
1 c NLSL Version 1.5 7/2/95
2 c---------------------------------------------------------------------
3 c =========================
4 c subroutine ORDRPR
5 c =========================
6 c
7 c This routine calculates the molecular order parameters corresponding
8 c to a specified axial orienting potential. The potential is
9 c expressed as
10 c
11 c L L
12 c U(phi,theta,0) = Sum Sum c D (phi,theta,0)
13 c L K K 0K
14 c
15 c where D^L_{0K} is a Wigner rotation matrix element (equivalent to
16 c the spherical harmonic Y^L_K), c^L_K are user-specified coefficients,
17 c and the sums are restricted to L=0 to 4 and K=0 to L (for even L,K only)
18 c
19 c The coefficients are specified in the cf array in the order
20 c c20, c22, c40, c42, c44. The routine calculates an order parameter for
21 c each nonzero coefficient specififed. The order parameter is the average
22 c of <D^L_{0K}(phi,theta)> weighted by the value of exp(-U(phi,theta))
23 c and integrated over angles theta and phi.
24 c
25 c
26 c Includes:
27 c rndoff.inc
28 c
29 c Uses:
30 c ccrints.f (Integration routines ccrint ccrin1)
31 c ftht (Function returning theta-dependent part of U)
32 c fphi (Function returning phi-dependent part of U)
33 c
34 c The last two functions are included in this file.
35 c Written by David E. Budil
36 c Cornell Department of Chemistry, March 1992
37 c
38 c----------------------------------------------------------------------
39 
40  subroutine ordrpr(cf,d)
41 c
42  use rnddbl
43 c
44  implicit none
45  double precision cf(5),d(5)
46 c
47  double precision ZERO,ONE,SMALL
48  parameter(zero=0.0d0,one=1.0d0,small=1.0d-16)
49 c
50  integer i,npt,nup,id
51  double precision fz,pfn
52 c
53  integer kount
54  double precision acc,c,d20,d40,bd22,bd42,bd44
55  common/potprm/acc,c(5),kount
56  common/func/d20,d40,bd22,bd42,bd44
57 c
58  double precision ftht
59  external ftht
60 c......................................................................
61  acc=0.001
62  npt=0
63 c
64  do kount=1,5
65  c(kount)=cf(kount)
66  if (dabs(c(kount)).gt.rndoff) npt=kount
67  end do
68 c
69  if (npt.gt.0) then
70 c
71 c--------------------------------------------------
72 c Integrate unweighted potential function (kount=0)
73 c--------------------------------------------------
74  kount=0
75  call ccrint(zero,one,acc,small,pfn,nup,ftht,id)
76 c
77 c-------------------------------------------------------------------
78 c Integrate potential weighted by D20,D22,D40,D42,D44 for kount=1..5
79 c Only calculate up to order of highest-order potential coefficient (npt)
80 c-------------------------------------------------------------------
81  do kount=1,npt
82  call ccrint(zero,one,acc,small,fz,nup,ftht,id)
83  d(kount)=fz/pfn
84  end do
85  endif
86  return
87  end
88 
89 c----------------------------------------------------------------------
90 c =========================
91 c function FTHT
92 c =========================
93 c
94 c Contains theta-dependence of orienting pseudopotential
95 c----------------------------------------------------------------------
96 c
97  function ftht(ctht)
98 c
99  use pidef
100 c
101  implicit none
102  double precision ftht,ctht
103 c
104  integer nup,id
105  double precision ctht2,stht2,result
106 c
107  integer kount
108  double precision acc,c,d20,d40,bd22,bd42,bd44
109  common/potprm/acc,c(5),kount
110  common/func/d20,d40,bd22,bd42,bd44
111 c
112  double precision fphi
113  external fphi
114 c
115 c---------------------------------------------------------
116 c definition of some constants
117 c A22 = sqrt(3/2), A42 = sqrt(5/8), A44 = sqrt(35/8)
118 c---------------------------------------------------------
119  double precision A22,A42,A44,ONE,ZERO,SMALL
120  parameter(a22=1.22474487139159d0,
121  2 a42=0.790569415042095d0,
122  3 a44=1.04582503316759d0 )
123  parameter(one=1.0d0,zero=0.0d0,small=1.0d-16 )
124 c
125 c......................................................................
126 c
127  ctht2=ctht*ctht
128  stht2=one-ctht2
129  d20 =1.5d0*ctht2-0.5d0
130  d40 =( (4.375d0*ctht2)-3.75d0)*ctht2+0.375d0
131  bd22=a22*stht2
132  bd42=a42*stht2*(7.0d0*ctht2-one)
133  bd44=a44*stht2*stht2
134 c
135  call ccrin1(zero,pi,acc,small,result,nup,fphi,id)
136  ftht=result
137  return
138  end
139 
140 c----------------------------------------------------------------------
141 c =========================
142 c function FPHI
143 c =========================
144 c
145 c Contains phi-dependence of orienting pseudopotential
146 c----------------------------------------------------------------------
147  function fphi(phi)
148  implicit none
149  double precision fphi,phi,c2phi,c4phi
150 c
151  integer kount
152  double precision acc,c,d20,d40,bd22,bd42,bd44
153  common/potprm/acc,c(5),kount
154  common/func/d20,d40,bd22,bd42,bd44
155 c
156  double precision ONE,TWO
157  parameter( one=1.0d0,
158  # two=2.0d0 )
159 c
160  c2phi=dcos(phi+phi)
161  c4phi=two*c2phi*c2phi - one
162 c
163  fphi=dexp( c(1)*d20
164  2 + c(2)*bd22*c2phi
165  3 + c(3)*d40
166  4 + c(4)*bd42*c2phi
167  5 + c(5)*bd44*c4phi )
168 c
169  if(kount.eq.0) return
170  if(kount.eq.1) fphi=d20*fphi
171  if(kount.eq.2) fphi=bd22*fphi*c2phi
172  if(kount.eq.3) fphi=d40*fphi
173  if(kount.eq.4) fphi=bd42*fphi*c2phi
174  if(kount.eq.5) fphi=bd44*fphi*c4phi
175  return
176  end
177 
178 c----------------------------------------------------------------------
179 c =========================
180 c function FU20
181 c =========================
182 c Function defining Boltzmann population at a given orientation
183 c resulting from an orienting pseudopotential using only the c20
184 c term. Used for calculating weighting as a function of angle
185 c between the average director and the spectrometer field, for
186 c the "partial MOMD" model. The potential is intended to describe
187 c partial macroscopic orientation of the directors of local domains.
188 c
189 c This is required to be a single-argument function so that it
190 c may be called by ccrint to calculate the normalization factor
191 c for orientational weighting coefficients. The second parameter, the
192 c director c20 coefficient, is passed in common /dfunc/
193 c----------------------------------------------------------------------
194  function fu20(cost)
195 c
196  use dfunc
197 c
198  implicit none
199  double precision fu20,cost
200 c
201  fu20 = exp( dlam*(1.5d0*cost*cost-0.5d0 ) )
202  return
203  end
204 c
205 c
206 c phi is expected in radians
207 c
208  function fu20phi(phi)
209 c
210  use dfunc
211  use pidef
212 c
213  implicit none
214  double precision costht,fu20phi,phi
215 c
216  costht = ss*cos(phi)+cc
217  fu20phi = exp( dlam*(1.5d0*costht*costht-0.5d0 ) )
218  return
219  end
double precision function fu20phi(phi)
Definition: ordrpr.f90:209
double precision function fphi(phi)
Definition: ordrpr.f90:148
subroutine ccrin1(bndlow, bndhi, epsiln, small, sum, neval, f, id)
Definition: ccrints.f90:240
double precision, save cc
Definition: dfunc.f90:19
double precision, save ss
Definition: dfunc.f90:19
double precision function fu20(cost)
Definition: ordrpr.f90:195
double precision function ftht(ctht)
Definition: ordrpr.f90:98
subroutine ordrpr(cf, d)
Definition: ordrpr.f90:41
subroutine ccrint(bndlow, bndhi, epsiln, small, sum, neval, f, id)
Definition: ccrints.f90:32
double precision, parameter pi
Definition: pidef.f90:15
double precision, save dlam
Definition: dfunc.f90:19
Definition: dfunc.f90:16
double precision, parameter rndoff
Definition: rnddbl.f90:86
Definition: pidef.f90:12