NLSL
matrll.f90
Go to the documentation of this file.
1 c NLSL Version 1.5.1 beta 11/25/95
2 c**********************************************************************
3 c =========================
4 c subroutine MATRLL
5 c =========================
6 c
7 c Subroutine for calculating matrix elements for EPRLL calculation.
8 c
9 c Updated to improve matrix storage algorithm. This version differs
10 c from the published version (Biol. Magnetic Resonance, Vol. 8)
11 c in that it only stores half of the matrix. Changes are noted
12 c in the code given below.
13 c
14 c Updated to provide for a nonaxial diffusion tensor R
15 c (only in the case of Brownian motion).
16 c
17 c Updated to replace input of residence times tl, tkxy, tkzz
18 c with the product of residence time and diffusion rate constants:
19 c pml, pmzz, pmxy. (Differs from parameters used in the EPRLL/EPRCGL
20 c programs in this respect)
21 c
22 c flags
23 c -----
24 c
25 c fld : True if |lr-lc| <=2 , false otherwise.
26 c Used as flag for bandwidth of Liouville operator
27 c matrix elements. The restriction on the
28 c L-bandwidth due the inclusion of spin operators
29 c of ranks 0 and 2.
30 c
31 c fldmd : True if |lr-lc| <=2 or mr=mc, false otherwise.
32 c Used as flag for possible non-zero matrix
33 c elements. The L-bandwidth restriction is
34 c explained above and the diffusion operator
35 c matrix elements can be found only when mr=mc.
36 c
37 c flk0 : True if lr=lc and kr=kc, false otherwise.
38 c Used as a flag for determing the existence of
39 c non-zero Liouville (A-tensor) and diffusion
40 c (Heisenberg spin exchange) matrix elements.
41 c
42 c written by DJS 29-AUG-87
43 c
44 c modified by JPG 06-JUN-88 to correct the anisotropic
45 c viscosity terms in the diffusion
46 c operator to conform to Polnaszek
47 c and Freed 1975 paper
48 c modified by DJS JUN-90 to update matrix storage. Comments by DEB.
49 c
50 c modified by DEB: calculation of isotropic part of diffusion
51 c operator moved outside iqnr loop. Fully anisotropic
52 c rotational diffusion tensor implemented.
53 c
54 c Uses:
55 c cd2km.f
56 c anxlk.f
57 c ipar.f
58 c w3j.f
59 c
60 c**********************************************************************
61 c
62  subroutine matrll(ierr)
63 c
64  use nlsdim
65  use rnddbl
66  use eprmat
67  use eprprm
68  use errmsg
69  use physcn
70  use stdio
71  use maxl
72 c
73  implicit none
74  integer i,j,l,kdi,ksi,li,ierr
75  integer lr,lrprod,iparlr,kr,krmx,mr,mrmx,
76  # ipnr,ipnrmn,ipnrmx,iqnr,iqnrmn,iqnrmx,
77  # lc,lcmax,iparlc,ld,ldabs,kc,kcmn,kcmx,kd,ks,
78  # kdabs,ksabs,mc,mcmn,mcmx,md,ms,mdabs,msabs,
79  # ipnc,ipncmn,ipncmx,iqnc,iqncmn,iqncmx,
80  # ipnd,ipns,ipndab,ipnsab,iqnd,iqns,iqndab,
81  # kip,nrow,ncol,nelr,neli,nel,ipnrsg,ipncsg,
82  # jmc,jmr,krmn,mrmn,jkc,jkr,jkd,jmd,kcsgn,krsgn,mcsgn,mrsgn
83 c
84  double precision dsq2,dsq6,dsq15,dsq23
85  double precision ct,rpl,rmi,rz,rj,rp,ossc,fii,ra0,
86  # cpmkr,cnpr,c1,c2,c3,c4,
87  # cdfd,cnl,cplkc,cnk,z,ra1,rg1,rw1,ra2,rg2,rw2,
88  # ra,rg,rw,cliou,cgam,cplmc,wliou1,wliou2,cdff,cdffr,
89  # cdffu,cnpc,cnp,cnorm,d1,d2,sa1,sa2,sw1,sw2,clioi1,fki,
90  # clioi2,clioua,sg1,sg2,clioug,ctemp,cgamw
91  double precision cjkr,cjkc,cjmc,zeen
92  double precision d2km(2,5,5)
93 c
94  logical flk0,fld,fldmd,fdpqi,fdjkm
95 c
96  double precision ZERO,ONE
97  parameter(zero=0.0d0,one=1.0d0)
98 c
99  integer ipar
100  double precision w3j
101  external ipar,w3j
102 c
103 c######################################################################
104 c
105 c.... Debugging purposes only!
106 c open (unit=20,file='nlmat.tst',status='unknown',
107 c # access='sequential',form='formatted')
108 c............................
109 c
110 c----------------------------------------------------------------------
111 c Scale diffusion constants
112 c----------------------------------------------------------------------
113 c
114  ct=g0*betae/hbar
115  rpl=0.5d0*(dx+dy)/ct
116  rmi=0.25*(dx-dy)/ct
117  rz=dz/ct
118  rj=djf/ct
119  rp=djfprp/ct
120  ossc=oss/ct
121 c
122 c----------------------------------------------------------------------
123 c Define constants
124 c----------------------------------------------------------------------
125 c
126  dsq2=dsqrt(2.0d0)
127  dsq6=dsqrt(6.0d0)
128  dsq15=dsqrt(1.5d0)
129  dsq23=dsqrt(2.0d0/3.0d0)
130 c
131  fii=0.25d0*dble(in2*(in2+2))
132 c
133  if(ipdf.ne.1) then
134  pml=zero
135  pmzz=zero
136  pmxy=zero
137  end if
138 c
139  if ((ipdf.eq.2) .and. (ist.eq.0)) then
140  rpl=rpl+rp
141  rz=rz+rp
142  end if
143 c
144  if(in2.eq.0) then
145  ra0=zero
146  zeen=zero
147  else
148  ra0=-dsq15*a0
149  zeen=b0*gamman/gammae
150  end if
151 c
152 c----------------------------------------------------------------------
153 c 2
154 c Calculate d (psi) as in Brink and Satchler
155 c k,m
156 c 2
157 c Note : dlm(i+3,j+3)=d (psi)
158 c i,j
159 c----------------------------------------------------------------------
160 c
161  call cd2km(d2km,zero,psi,zero)
162 c
163 c----------------------------------------------------------------------
164 c L
165 c Calculate the quantities X used in potential-dependent
166 c K
167 c portion of diffusion operator
168 c----------------------------------------------------------------------
169 c
170  call anxlk(rpl,rmi,rz)
171 c
172 c----------------------------------------------------------------------
173 c Initialize counters
174 c----------------------------------------------------------------------
175 c
176  nrow=0
177  ncol=0
178  nelr=0
179  neli=0
180  nel=0
181 c
182 c----------------------------------------------------------------------
183 c *** loop over lr ***
184 c----------------------------------------------------------------------
185 c
186  do 300 lr=0,lemx,ldelta
187 c
188  lrprod=lr*(lr+1)
189  iparlr=ipar(lr)
190 c
191 c --- Skip odd L greater than odd L maximum
192 c
193  if((lr.gt.lomx).and.(iparlr.ne.1)) go to 300
194 c
195 c----------------------------------------------------------------------
196 c *** loop over kr ***
197 c----------------------------------------------------------------------
198 c
199 c --- Find limits for K index for row
200 c
201  krmx = min0(lr,kmx)
202  krmn = max0(-lr,kmn)
203  do 320 krsgn=krmn,krmx,kdelta
204 c
205 c --- Find symmetry index jK from sign of K
206 c
207  kr=abs(krsgn)
208  jkr=isign(1,krsgn)
209  if (krsgn.eq.0) jkr=iparlr
210  if (jkr.lt.jkmn) goto 320
211 c
212 c --------------
213 c Find parities
214 c --------------
215  cjkr=jkr
216 c
217 c
218 c----------------------------------------------------------------------
219 c *** loop over mr ***
220 c----------------------------------------------------------------------
221 c
222 c --- Find limits for M index for row
223 c
224  mrmx=min0(lr,mmx)
225  mrmn=max0(-lr,mmn)
226  do 340 mrsgn=mrmn,mrmx
227 c
228 c --- Find symmetry index jM from sign of M
229 c (if M=0, jM is specified by sign of pI)
230 c
231  mr=abs(mrsgn)
232  jmr=isign(1,mrsgn)
233 c
234 c --------------
235 c find parities
236 c --------------
237  cpmkr=ipar(mr+kr)
238 c
239 c----------------------------------------------------------------------
240 c Calculate isotropic part of diffusion operator
241 c (Code moved to here since it does not depend on nuclear
242 c transition indices--DEB)
243 c----------------------------------------------------------------------
244 c
245  c1=dble(lrprod)
246  if(ipdf.ne.0) then
247  c2=one+pml*c1
248  if(dabs(expl).gt.rndoff) c2=c2**expl
249  c1=c1/c2
250  end if
251 c
252  cdfd=rpl*c1
253  c1=dble(kr*kr)
254  c2=rz
255  c3=rpl
256 c
257  if(ipdf.ne.0) then
258  c4=one+pmzz*c1
259  if(dabs(expkzz).gt.rndoff) c4=c4**expkzz
260  c2=c2/c4
261  c4=one+pmxy*c1
262  if(dabs(expkxy).gt.rndoff) c4=c4**expkxy
263  c3=c3/c4
264  end if
265 c
266  cdfd=cdfd+c1*(c2-c3)
267  if(ipdf.eq.2) cdfd=cdfd+dble(mr*mr)*(rj-rp)
268 c
269  if(ist.ne.0) then
270  i=kr/ist
271  if((i*ist).ne.kr) cdfd=cdfd+rj
272  end if
273 c
274 c --------------------------------
275 c Find limits for pI index for row
276 c --------------------------------
277  if(ipsi0.eq.0) then
278  ipnrmn=mr
279  ipnrmx=mr
280  else
281  ipnrmx = min(in2,ipnmx)
282  if(mr.eq.0 .and. jmmn.eq.1) then
283  ipnrmn=0
284  else
285  ipnrmn=-ipnrmx
286  end if
287  end if
288 c
289 c----------------------------------------------------------------------
290 c *** loop over ipnr ***
291 c----------------------------------------------------------------------
292 c
293  do 350 ipnrsg=ipnrmn,ipnrmx
294 c
295 c --- Find symmetry index jM from sign of pI if undetermined
296 c
297  if (mr.eq.0) then
298  jmr=isign(1,ipnrsg)
299  ipnr=abs(ipnrsg)
300  if (ipnr.eq.0) jmr=iparlr
301  else
302  ipnr=ipnrsg
303  end if
304  if(jmr.lt.jmmn) goto 350
305 c
306 c----------------------------------------------------------------------
307 c Find scaling factor for matrix elements
308 c----------------------------------------------------------------------
309 c
310  if((ipnr.eq.0).and.(mr.eq.0)) then
311  cnpr=dsq2
312  else
313  cnpr=one
314  end if
315 c
316 c----------------------------------------------------------------------
317 c Find limits for qI index for row
318 c----------------------------------------------------------------------
319 c
320  iqnrmx=in2-abs(ipnr)
321  iqnrmn=-iqnrmx
322 c
323 c----------------------------------------------------------------------
324 c *** loop over iqnr ***
325 c----------------------------------------------------------------------
326 c
327  do 360 iqnr=iqnrmn,iqnrmx,2
328 c
329 c----------------------------------------------------------------------
330 c Increment row counter
331 c----------------------------------------------------------------------
332 c
333  nrow=nrow+1
334  if (nrow.gt.mxdim) then
335  ierr=dimbig
336  go to 9999
337  end if
338 c
339  jzmat(nrow)=neli+1
340  kzmat(nrow)=nelr+1
341 c
342 c----------------------------------------------------------------------
343 c Calculate limits for L index for columns
344 c----------------------------------------------------------------------
345 c
346  lcmax = min0( lr+lband, lemx )
347 c
348 c----------------------------------------------------------------------
349 c Initialize column counter
350 c----------------------------------------------------------------------
351 c
352  ncol=nrow
353 c
354 c Check for user input or halt
355 c
356  call wpoll
357  if (hltcmd.ne.0 .or. hltfit.ne.0) then
358  ierr=mtxhlt
359  return
360  end if
361 c
362 c----------------------------------------------------------------------
363 c *** loop over lc ***
364 c----------------------------------------------------------------------
365 c
366  do 400 lc=lr,lcmax,ldelta
367 c
368  iparlc=ipar(lc)
369 c
370 c----------------------------------------------------------------------
371 c skip if L is odd and greater than odd L maximum
372 c----------------------------------------------------------------------
373 c
374  if((lc.gt.lomx).and.(iparlc.ne.1)) go to 400
375 c
376 c----------------------------------------------------------------------
377 c find limits for non-zero Liouville operator
378 c matrix elements
379 c----------------------------------------------------------------------
380 c
381  ld=lr-lc
382  ldabs=abs(ld)
383  fld=abs(ld).le.2
384 c
385 c --- find scaling factor
386 c
387  cnl=dsqrt((2.d0*lr+1.d0)*(2.d0*lc+1.d0))
388 c
389 c --- Find limits of K index for columns
390 
391  kcmx=min(lc,kmx)
392 c
393  if (nrow.eq.ncol) then
394  kcmn=kr
395  else
396  kcmn=max(-lc,kmn)
397  endif
398 c
399 c----------------------------------------------------------------------
400 c *** loop over kc ***
401 c----------------------------------------------------------------------
402 c
403  do 420 kcsgn=kcmn,kcmx,kdelta
404 c
405 c --- Find symmetry index of jK from sign of K
406 c
407  kc=abs(kcsgn)
408  jkc=isign(1,kcsgn)
409  if (kcsgn.eq.0) jkc=iparlc
410  if(jkc.lt.jkmn) go to 420
411 c
412 c --- find sums, differences and parities
413 c
414  kd=kr-kc
415  ks=kr+kc
416  jkd=jkr-jkc
417  kdabs=abs(kd)
418  ksabs=abs(ks)
419  cplkc=ipar(lc+kc)
420  cjkc=jkc
421 c
422 c -- Find scaling factor
423 c
424  if((kr.eq.0).and.(kc.eq.0)) then
425  cnk=0.5d0
426  else if((kr.ne.0).and.(kc.ne.0)) then
427  cnk=one
428  else
429  cnk=one/dsq2
430  end if
431 c
432 c----------------------------------------------------------------------
433 c Calculate rhombic part of diffusion tensor
434 c----------------------------------------------------------------------
435 c
436  if (kr+2.eq.kc .and. jkd.eq.0) then
437  cdffr=rmi*sqrt(dble((lr+kc-1)*(lr+kc)*
438  # (lr-kc+1)*(lr-kc+2) ) )/cnk
439  else if (kr-2.eq.kc .and. jkd.eq.0) then
440  cdffr=rmi*sqrt(dble((lr-kc-1)*(lr-kc)*
441  # (lr+kc+1)*(lr+kc+2) ) )/cnk
442  else
443  cdffr=zero
444  end if
445 c
446 c---------------------------------------------------------------------
447 c Calculate R(mu,l) as in Meirovitch, Igner, Igner, Moro, & Freed
448 c J. Phys. Chem. 77 (1982) p. 3935, Eqs. A42 & A44
449 c mu = g, A, W
450 c---------------------------------------------------------------------
451 c
452  if((kdabs.le.2).and.fld) then
453  z=w3j(lr,2,lc,kr,-kd,-kc)
454  if (jkd.eq.0) then
455  ra1=fad(1,kd+3)*z
456  rg1=fgd(1,kd+3)*z
457  rw1=fwd(1,kd+3)*z
458  else
459  ra1=fad(2,kd+3)*z*cjkr
460  rg1=fgd(2,kd+3)*z*cjkr
461  rw1=fwd(2,kd+3)*z*cjkr
462  end if
463  else
464  ra1=zero
465  rg1=zero
466  rw1=zero
467  end if
468 c
469  if((ksabs.le.2).and.fld) then
470  if (jkd.eq.0) then
471  z=w3j(lr,2,lc,kr,-ks,kc)
472  ra2=fad(1,ks+3)*z
473  rg2=fgd(1,ks+3)*z
474  rw2=fwd(1,ks+3)*z
475  else
476  ra2=fad(2,ks+3)*z*cjkr
477  rg2=fgd(2,ks+3)*z*cjkr
478  rw2=fwd(2,ks+3)*z*cjkr
479  end if
480  else
481  ra2=zero
482  rg2=zero
483  rw2=zero
484  end if
485 c
486  if(in2.ne.0) then
487  ra=ra1+cplkc*cjkc*ra2
488  else
489  ra=zero
490  end if
491 c
492  rg=rg1+cplkc*cjkc*rg2
493  rw=rw1+cplkc*cjkc*rw2
494 c
495 c----------------------------------------------------------------------
496 c find limits for M index for columns
497 c----------------------------------------------------------------------
498 c
499  mcmx=min(lc,mmx)
500  if (nrow.eq.ncol) then
501  mcmn=mr
502  else
503  mcmn=max(-lc,mmn)
504  endif
505 c
506 c----------------------------------------------------------------------
507 c *** loop over mc ***
508 c----------------------------------------------------------------------
509 c
510  do 440 mcsgn=mcmn,mcmx
511 c
512 c --- Find symmetry index jM from sign of M
513 c (or pI if M is 0)
514 c
515  mc=abs(mcsgn)
516  jmc=isign(1,mcsgn)
517 c
518 c --- Find M differences, sums, and parities
519 c
520  md=mr-mc
521  ms=mr+mc
522  mdabs=abs(md)
523  msabs=abs(ms)
524  cplmc=ipar(lc+mc)
525 c
526 c --------------------------------------
527 c set flags and zero out matrix elements
528 c if no non-zero elements are possible
529 c --------------------------------------
530  flk0 = (ld.eq.0).and.(kd.eq.0)
531  if(fld.or.(md.eq.0)) then
532  fldmd=.true.
533  else
534  cliou=zero
535  cgam=zero
536  fldmd=.false.
537  end if
538 c
539 c -------------------------------------------------
540 c get appropriate 3-J symbols if non-zero Liouville
541 c operator matrix elements are possible
542 c -------------------------------------------------
543  if(fld) then
544  wliou1=w3j(lr,2,lc,mr,-md,-mc)
545  wliou2=w3j(lr,2,lc,mr,-ms,mc)
546  else
547  wliou1=zero
548  wliou2=zero
549  end if
550 c
551 c----------------------------------------------------------------------
552 c Calculate off-diagonal terms of the diffusion operator, including
553 c the potential-dependent portion and the rhombic component of the
554 c isotropic diffusion operator.
555 c
556 c See Meirovitch,Igner,Igner,Moro & Freed, J. Chem. Phys. 77 (1982)
557 c pp.3933-3934 especially Eqns. A22-24 and A40 for more information
558 c----------------------------------------------------------------------
559 c
560 c --- Rhombic part of isotropic diffusion operator
561 c
562  if (ld.eq.0 .and. md.eq.0) then
563  cdff=cdffr
564  else
565  cdff=zero
566  end if
567 c
568 c --- Potential-dependent part of diffusion operator
569 c
570  if((ipt.ne.0).and.(ldabs.le.lband).and.(md.eq.0)
571  # .and.((kdabs.le.kband).or.(ksabs.le.kband)).and.
572  # (ipar(ks).eq.1).and.(jkd.eq.0)) then
573 c
574  kdi=1+kdabs/2
575  ksi=1+ksabs/2
576  c1=cpmkr*cnl*cnk
577 c
578  do l=0,lband,2
579  cdffu=zero
580  li=1+l/2
581  if (ksi.le.li .and. abs(xlk(li,ksi)).ge.rndoff)
582  # cdffu=cjkc*cplkc*xlk(li,ksi)
583  # *w3j(lr,l,lc,kr,-ks,kc)
584 c
585  if (kdi.le.li .and. abs(xlk(li,kdi)).ge.rndoff)
586  # cdffu=cdffu+xlk(li,kdi)*w3j(lr,l,lc,kr,-kd,-kc)
587 c
588  if (abs(cdffu).gt.rndoff) cdff=cdff+
589  # w3j(lr,l,lc,mr,0,-mr)*c1*cdffu
590  end do
591 c
592  end if
593 c
594 c----------------------------------------------------------------------
595 c Find limits for pI index for columns
596 c----------------------------------------------------------------------
597 c
598  if(ipsi0.eq.0) then
599  ipncmn=mc
600  ipncmx=mc
601  else
602  ipncmx=min(ipnmx,in2)
603 c
604  if (mc.eq.0.and.jmmn.eq.1) then
605  ipncmn=0
606  else
607  ipncmn=-ipncmx
608  end if
609  end if
610 c
611  if (nrow.eq.ncol) ipncmn=ipnr
612 c
613 c----------------------------------------------------------------------
614 c *** loop over ipnc ***
615 c----------------------------------------------------------------------
616 c
617  do 450 ipncsg=ipncmn,ipncmx
618 c
619  if(mc.eq.0) then
620  ipnc=abs(ipncsg)
621  jmc=isign(1,ipncsg)
622  if (ipnc.eq.0) jmc=iparlc
623  else
624  ipnc=ipncsg
625  end if
626  if (jmc.lt.jmmn) goto 450
627 c
628 c --- Find pI sums and differences
629 c
630  ipnd=ipnr-ipnc
631  ipns=ipnr+ipnc
632  ipndab=abs(ipnd)
633  ipnsab=abs(ipns)
634 c
635 c ---find scaling factors
636 c
637  if((ipnc.eq.0).and.(mc.eq.0)) then
638  cnpc=dsq2
639  else
640  cnpc=one
641  end if
642 c
643  cjmc=jmc
644  jmd=jmr-jmc
645  fdjkm=(jkd.eq.0).and.(jmd.eq.0)
646  cnp=one/(cnpr*cnpc)
647  cnorm=cnl*cnk*cnp*cpmkr
648 c
649 c----------------------------------------------------------------------
650 c get director tilt dependence of Liouville
651 c operator matrix elements
652 c----------------------------------------------------------------------
653 c
654  if((ipndab.le.2).and.(mdabs.le.2)) then
655  d1=d2km(1,ipnd+3,md+3)*wliou1
656  else
657  d1=zero
658  end if
659 c
660  if((ipnsab.le.2).and.(msabs.le.2)) then
661  d2=d2km(1,ipns+3,ms+3)*wliou2
662  else
663  d2=zero
664  end if
665 c
666 c----------------------------------------------------------------------
667 c find limits to qI index for columns
668 c----------------------------------------------------------------------
669 c
670  iqncmx=in2-abs(ipnc)
671  iqncmn=-iqncmx
672 c
673 c
674 c----------------------------------------------------------------------
675 c *** loop over iqnc ***
676 c----------------------------------------------------------------------
677 c
678  if (nrow.eq.ncol) iqncmn=iqnr
679  do 460 iqnc=iqncmn,iqncmx,2
680 c
681 c----------------------------------------------------------------------
682 c skip calculation of matrix elements if non-zero
683 c matrix elements are not possible
684 c----------------------------------------------------------------------
685 c
686  if(.not.fldmd) go to 470
687 c
688 c----------------------------------------------------------------------
689 c find sums and differences
690 c----------------------------------------------------------------------
691 c
692  iqnd=iqnr-iqnc
693  iqns=iqnr+iqnc
694  iqndab=abs(iqnd)
695  fdpqi=(ipnd.eq.0).and.(iqnd.eq.0)
696 c
697 c----------------------------------------------------------------------
698 c Liouville operator
699 c
700 c NOTE: see Appendix B p.3936-3937 in M,I,I,M, & F for
701 c details. Sa and Sg in program also contain the
702 c appropriate Clebsch-Gordon coefficient. fki is
703 c Ki in text.
704 c The isotropic contributions from the A tensor are
705 c not described in the text. They appear only when
706 c L1 (lr) equals L2 (lc). The Clebsch-Gordon coefficients
707 c for these terms (+/- 1/sqrt(3)) are included in the
708 c definition of a0. The Wigner 3J symbols appropriate for
709 c the isotropic terms exactly cancel the normalization
710 c constant cnorm. Instead of calculating the 3J symbols,
711 c the cnorm factor is simply omitted from the isotropic
712 c terms below.
713 c
714 c----------------------------------------------------------------------
715 c
716  if(fld) then
717 c
718  if (jmd.eq.0) then
719 c
720  if((ipndab.eq.iqndab).and.
721  # (ipndab.le.1).and.(in2.ne.0)) then
722 c
723 c ----------------------------
724 c Hyperfine interaction term
725 c ----------------------------
726  if(ipnd.eq.0) then
727  sa1=iqnr/dsq6
728  if(flk0.and.md.eq.0.and.jkd.eq.0) then
729  clioi1=-sa1*ra0
730  else
731  clioi1=zero
732  end if
733 c
734  else
735  kip=iqnr*iqnd+ipnr*ipnd
736  kip=kip*(kip-2)
737  fki=dsqrt(fii-0.25d0*kip)
738  sa1=-ipnd*fki*0.25d0
739  clioi1=zero
740  end if
741 c
742  else
743  sa1=zero
744  clioi1=zero
745  end if
746 c
747  if((ipnsab.eq.iqndab).and.
748  # (ipnsab.le.1).and.(in2.ne.0)) then
749 c
750  if(ipns.eq.0) then
751  sa2=iqnr/dsq6
752  if(flk0.and.ms.eq.0.and.jkd.eq.0) then
753  clioi2=-sa2*ra0
754  else
755  clioi2=zero
756  end if
757  else
758  kip=iqnr*iqnd+ipnr*ipns
759  kip=kip*(kip-2)
760  fki=dsqrt(fii-0.25d0*kip)
761  sa2=-ipns*fki*0.25d0
762  clioi2=zero
763  end if
764 c
765  else
766  sa2=zero
767  clioi2=zero
768  end if
769 c
770  clioua=(sa1*d1+cjmc*cplmc*sa2*d2)
771 c
772 c -----------------------
773 c Electronic Zeeman term
774 c -----------------------
775  if((iqnd.eq.0).and.
776  # (abs(rg).gt.rndoff)) then
777 c
778  if(ipnd.eq.0) then
779  sg1=dsq23
780  else
781  sg1=zero
782  end if
783 c
784  if(ipns.eq.0) then
785  sg2=dsq23
786  else
787  sg2=zero
788  end if
789 c
790  clioug=(sg1*d1+cjmc*cplmc*sg2*d2)
791 c
792  else
793  clioug=zero
794  end if
795 c
796 c --------------------------------------
797 c Orientation-dependent linewidth term
798 c --------------------------------------
799  if((iqnd.eq.0).and.(abs(rw).gt.rndoff))
800  # then
801  if(ipnd.eq.0) then
802  sw1=dsq23
803  else
804  sw1=zero
805  end if
806 c
807  if(ipns.eq.0) then
808  sw2=dsq23
809  else
810  sw2=zero
811  end if
812  cgamw=cnorm*(sw1*d1+cjmc*cplmc*sw2*d2)
813  # *rw
814 c
815  else
816  cgamw=zero
817  end if
818 c
819  cliou=cnorm*(clioua*ra+clioug*rg)+
820  # cnp*(clioi1+clioi2)
821 c
822 c if (jmd.eq.0) then...
823  else
824 c
825 c --------------------
826 c Nuclear Zeeman term
827 c --------------------
828  if(flk0.and.(md.eq.0).and.(jkd.eq.0)
829  # .and.fdpqi) then
830  cliou=ipnr*zeen
831  else
832  cliou=zero
833  end if
834 c
835  cgamw=zero
836  end if
837 c
838 c if (fld) then...
839  else
840  cliou=zero
841  cgamw=zero
842  end if
843 c
844 c----------------------------------------------------------------------
845 c Heisenberg exchange operator
846 c Now also includes orientational averaging via HE
847 c in slow-motional region with the delta(L1,0) term.
848 c (See Lee,Budil,and Freed, JCP, 1994, Eq. B9)
849 c----------------------------------------------------------------------
850 c
851  if((abs(ossc).gt.rndoff).and.(ipnd.eq.0).and.
852  # (md.eq.0).and.flk0.and.fdjkm) then
853 c
854  ctemp=zero
855  if(iqnd.eq.0) ctemp=one
856  if ((ipnr.eq.0) .and.(lr.eq.0))
857  # ctemp=ctemp-one/dble(in2+1)
858  cgam=cgamw+ctemp*ossc
859 c
860  else
861  cgam=cgamw
862  end if
863 c
864 c ----------------------------------------
865 c Add in diffusion terms on the diagonal
866 c ----------------------------------------
867  if(fdpqi.and.fdjkm) cgam=cgam+cdff
868 c
869 c----------------------------------------------------------------------
870 c Store matrix elements
871 c Real off-diagonal matrix elements (cgam nonzero) are stored
872 c sequentially starting from the end of the space allocated for zmat
873 c and working down.
874 c
875 c Imaginary off-diagonal matrix elements (cliou nonzero) are stored
876 c sequentially starting at the beginning of the zmat array.
877 c----------------------------------------------------------------------
878 c
879 c --- Diagonal elements
880 c
881  if (nrow.eq.ncol) then
882  cgam=cgam+cdfd
883  zdiag(1,nrow)=cgam
884  zdiag(2,nrow)=-cliou
885 c
886 c --- Off-diagonal elements
887 c
888  else
889  if (abs(cliou).gt.rndoff) then
890  nel=nel+1
891  if (nel.gt.mxel) then
892  ierr=mtxbig
893  go to 9999
894  else
895  neli=neli+1
896  zmat(neli)=-cliou
897  izmat(neli)=ncol
898  end if
899  end if
900 c
901  if (abs(cgam).gt.rndoff) then
902  nel=nel+1
903  if (nel.gt.mxel) then
904  ierr=mtxbig
905  go to 9999
906  else
907  nelr=nelr+1
908  zmat(mxel-nelr+1)=cgam
909  izmat(mxel-nelr+1)=ncol
910  end if
911  end if
912 c
913 c if (ncol.eq.nrow)...else...
914  end if
915 c
916 c..........Debugging purposes only!
917 c if (abs(cgam).gt.RNDOFF.or.abs(cliou).gt.RNDOFF)
918 c # write (20,7000) lr,krsgn,mrsgn,ipnr,iqnr,
919 c # lc,kcsgn,mcsgn,ipnc,iqnc,cgam,-cliou
920 c 7000 format('<',4(i2,','),i2,'|L|',4(i2,','),i2,'> =',2g14.7)
921 c...................................
922 c
923 c
924 c --------------------------
925 c increment column counter
926 c --------------------------
927 c
928  470 ncol=ncol+1
929 c
930  if (ncol.gt.mxdim) then
931  ierr=dimbig
932  go to 9999
933  end if
934 c
935 c----------------------------------------------------------------------
936 c end loop over columns
937 c----------------------------------------------------------------------
938 c
939  460 continue
940  450 continue
941  440 continue
942  420 continue
943  400 continue
944 c
945 c----------------------------------------------------------------------
946 c end loop over rows
947 c----------------------------------------------------------------------
948 c
949  360 continue
950  350 continue
951  340 continue
952  320 continue
953  300 continue
954 c
955  ierr=0
956  ndim=nrow
957  neltot=nel
958  nelre=nelr
959  nelim=neli
960 c
961  jzmat(ndim+1)=neli+1
962  kzmat(ndim+1)=nelr+1
963 c
964 c.......... Debugging purposes only!
965 c close (20)
966 c...................................
967  9999 return
968  end
double precision, pointer, save pml
Definition: eprprm.f90:58
integer, save ipt
Definition: eprprm.f90:82
double precision, pointer, save dz
Definition: eprprm.f90:58
double precision, pointer, save dx
Definition: eprprm.f90:58
void FORTRAN() wpoll()
Definition: pltx.c:860
subroutine cd2km(d2km, alpha, beta, gamma)
Definition: cd2km.f90:33
double precision, save expkzz
Definition: eprprm.f90:78
integer, dimension(mxdim+1), save jzmat
Definition: eprmat.f90:50
Definition: maxl.f90:22
double precision, save expl
Definition: eprprm.f90:78
integer, save lband
Definition: eprprm.f90:82
integer, pointer, save ist
Definition: eprprm.f90:69
integer, pointer, save ndim
Definition: eprprm.f90:69
integer, save hltcmd
Definition: stdio.f90:32
double precision, parameter betae
Definition: physcn.f90:20
double precision, dimension(2, 5), save fad
Definition: eprprm.f90:78
subroutine anxlk(rp, rm, rz)
Definition: anxlk.f90:37
integer, parameter mtxbig
Definition: errmsg.f90:51
integer, save hltfit
Definition: stdio.f90:32
Definition: stdio.f90:26
integer, pointer, save mmx
Definition: eprprm.f90:69
integer, pointer, save ipdf
Definition: eprprm.f90:69
double precision, pointer, save dy
Definition: eprprm.f90:58
integer, pointer, save jkmn
Definition: eprprm.f90:69
double precision, pointer, save gamman
Definition: eprprm.f90:58
double precision, dimension(2, 5), save fgd
Definition: eprprm.f90:78
integer, dimension(mxdim+1), save kzmat
Definition: eprmat.f90:50
integer, dimension(mxel), save izmat
Definition: eprmat.f90:50
integer, pointer, save mmn
Definition: eprprm.f90:69
integer, pointer, save ipnmx
Definition: eprprm.f90:69
integer, save nelre
Definition: eprprm.f90:82
integer, save nelim
Definition: eprprm.f90:82
double precision, pointer, save psi
Definition: eprprm.f90:58
double precision, pointer, save pmxy
Definition: eprprm.f90:58
integer, save kband
Definition: eprprm.f90:82
double precision, save expkxy
Definition: eprprm.f90:78
double precision, pointer, save djfprp
Definition: eprprm.f90:58
double precision, dimension(mxel), save zmat
Definition: eprmat.f90:49
integer, save ldelta
Definition: eprprm.f90:82
integer, pointer, save lomx
Definition: eprprm.f90:69
integer, pointer, save kmn
Definition: eprprm.f90:69
double precision, dimension(5, 5), save xlk
Definition: eprprm.f90:78
integer, parameter mxel
Definition: nlsdim.f90:39
double precision, pointer, save pmzz
Definition: eprprm.f90:58
integer, pointer, save in2
Definition: eprprm.f90:69
integer, pointer, save kmx
Definition: eprprm.f90:69
subroutine matrll(ierr)
Definition: matrll.f90:63
double precision, save g0
Definition: eprprm.f90:78
integer, pointer, save jmmn
Definition: eprprm.f90:69
integer, pointer, save lemx
Definition: eprprm.f90:69
double precision, parameter gammae
Definition: physcn.f90:20
integer, parameter mtxhlt
Definition: errmsg.f90:51
double precision, parameter hbar
Definition: physcn.f90:20
double precision, pointer, save b0
Definition: eprprm.f90:58
integer, parameter mxdim
Definition: nlsdim.f90:39
integer, parameter dimbig
Definition: errmsg.f90:51
integer, save ipsi0
Definition: eprprm.f90:82
double precision, dimension(2, 5), save fwd
Definition: eprprm.f90:78
double precision, pointer, save oss
Definition: eprprm.f90:58
double precision, save a0
Definition: eprprm.f90:78
integer, save neltot
Definition: eprprm.f90:82
double precision, parameter rndoff
Definition: rnddbl.f90:86
integer, save kdelta
Definition: eprprm.f90:82
double precision, dimension(2, mxdim), save zdiag
Definition: eprmat.f90:49
double precision, pointer, save djf
Definition: eprprm.f90:58