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