387 lines
9.1 KiB
Fortran
387 lines
9.1 KiB
Fortran
module ctrans_mod
|
|
use dim_parameter, only: qn
|
|
contains
|
|
!! subroutine ctrans
|
|
|
|
subroutine ctrans(q,x1,x2)
|
|
implicit none
|
|
include 'nnparams.incl'
|
|
include 'JTmod.incl'
|
|
double precision,intent(in):: q(qn)
|
|
double precision,intent(out):: x1(qn),x2(qn)
|
|
double precision:: cart(3,4),qint(maxnin)
|
|
integer i
|
|
!cart(:,1)=0.0d0
|
|
|
|
!cart(1:3,2:4) = reshape([ q(4:12) ], shape(cart(1:3,2:4)))
|
|
cart(1,1)=q(1)
|
|
cart(2,1)=q(2)
|
|
cart(3,1)=q(3)
|
|
cart(1,2)=q(4)
|
|
cart(2,2)=q(5)
|
|
cart(3,2)=q(6)
|
|
cart(1,3)=q(7)
|
|
cart(2,3)=q(8)
|
|
cart(3,3)=q(9)
|
|
cart(1,4)=q(10)
|
|
cart(2,4)=q(11)
|
|
cart(3,4)=q(12)
|
|
call cart2int(cart,qint)
|
|
do i=1,qn
|
|
if (abs(qint(i)) .lt. 1.0d-5) qint(i) =0.0d0
|
|
enddo
|
|
x1(1:qn)=qint(1:qn)
|
|
!x1(2)=0.0d0
|
|
x1(5)=-x1(5)
|
|
x1(3)=-x1(3)
|
|
!x1(6)=0.0d0
|
|
x2(1:qn)=0.0d0 !qint(1:qn)
|
|
end subroutine ctrans
|
|
|
|
|
|
subroutine cart2int(cart,qint)
|
|
implicit none
|
|
! This version merges both coordinate transformation routines into
|
|
! one. JTmod's sscales(2:3) are ignored.
|
|
! This is the first version to be compatible with one of my proper 6D fits
|
|
! Time-stamp: <2024-10-22 13:52:59 dwilliams>
|
|
|
|
! Input (cartesian, in Angström)
|
|
! cart(:,1): N
|
|
! cart(:,1+i): Hi
|
|
! Output
|
|
! qint(i): order defined in JTmod.
|
|
! Internal Variables
|
|
! no(1:3): NO distances 1-3
|
|
! pat_in: temporary coordinates
|
|
! axis: main axis of NO3
|
|
include 'nnparams.incl'
|
|
include 'JTmod.incl'
|
|
|
|
|
|
real*8 cart(3,4),qint(maxnin)
|
|
|
|
real*8 no(3), r1, r2, r3
|
|
real*8 v1(3), v2(3), v3(3)
|
|
real*8 n1(3), n2(3), n3(3), tr(3)
|
|
real*8 ortho(3)
|
|
real*8 pat_in(maxnin)
|
|
logical ignore_umbrella,dbg_umbrella
|
|
logical dbg_distances
|
|
|
|
!.. Debugging parameters
|
|
!.. set umbrella to 0
|
|
parameter (ignore_umbrella=.false.)
|
|
! parameter (ignore_umbrella=.true.)
|
|
|
|
!.. break if umbrella is not 0
|
|
parameter (dbg_umbrella=.false.)
|
|
! parameter (dbg_umbrella=.true.)
|
|
|
|
!.. break for tiny distances
|
|
parameter (dbg_distances=.false.)
|
|
! parameter (dbg_distances=.true.)
|
|
|
|
integer k
|
|
|
|
!.. get N-O vectors and distances:
|
|
do k=1,3
|
|
v1(k)=cart(k,2)-cart(k,1)
|
|
v2(k)=cart(k,3)-cart(k,1)
|
|
v3(k)=cart(k,4)-cart(k,1)
|
|
enddo
|
|
no(1)=norm(v1,3)
|
|
no(2)=norm(v2,3)
|
|
no(3)=norm(v3,3)
|
|
|
|
!.. temporarily store displacements
|
|
do k=1,3
|
|
pat_in(k)=no(k)-offsets(1)
|
|
enddo
|
|
|
|
do k=1,3
|
|
v1(k)=v1(k)/no(1)
|
|
v2(k)=v2(k)/no(2)
|
|
v3(k)=v3(k)/no(3)
|
|
enddo
|
|
|
|
!.. compute three normal vectors for the ONO planes:
|
|
call xprod(n1,v1,v2)
|
|
call xprod(n2,v2,v3)
|
|
call xprod(n3,v3,v1)
|
|
|
|
do k=1,3
|
|
tr(k)=(n1(k)+n2(k)+n3(k))/3.d0
|
|
enddo
|
|
r1=norm(tr,3)
|
|
do k=1,3
|
|
tr(k)=tr(k)/r1
|
|
enddo
|
|
|
|
! rotate trisector
|
|
call rot_trisec(tr,v1,v2,v3)
|
|
|
|
!.. determine trisector angle:
|
|
if (ignore_umbrella) then
|
|
pat_in(7)=0.0d0
|
|
else
|
|
pat_in(7)=pi/2.0d0 - acos(scalar(v1,tr,3))
|
|
pat_in(7)=sign(pat_in(7),cart(1,2))
|
|
endif
|
|
|
|
!.. molecule now lies in yz plane, compute projected ONO angles:
|
|
v1(1)=0.d0
|
|
v2(1)=0.d0
|
|
v3(1)=0.d0
|
|
r1=norm(v1,3)
|
|
r2=norm(v2,3)
|
|
r3=norm(v3,3)
|
|
do k=2,3
|
|
v1(k)=v1(k)/r1
|
|
v2(k)=v2(k)/r2
|
|
v3(k)=v3(k)/r3
|
|
enddo
|
|
|
|
! make orthogonal vector to v3
|
|
ortho(1)=0.0d0
|
|
ortho(2)=v3(3)
|
|
ortho(3)=-v3(2)
|
|
|
|
!.. projected ONO angles in radians
|
|
pat_in(4)=get_ang(v2,v3,ortho)
|
|
pat_in(5)=get_ang(v1,v3,ortho)
|
|
|
|
pat_in(6)=dabs(pat_in(5)-pat_in(4))
|
|
|
|
!.. account for rotational order of atoms
|
|
if (pat_in(4).le.pat_in(5)) then
|
|
pat_in(5)=2*pi-pat_in(4)-pat_in(6)
|
|
else
|
|
pat_in(4)=2*pi-pat_in(5)-pat_in(6)
|
|
endif
|
|
|
|
pat_in(4)=rad2deg*pat_in(4)-offsets(2)
|
|
pat_in(5)=rad2deg*pat_in(5)-offsets(2)
|
|
pat_in(6)=rad2deg*pat_in(6)-offsets(2)
|
|
pat_in(7)=rad2deg*pat_in(7)
|
|
|
|
call genANN_ctrans(pat_in)
|
|
|
|
qint(:)=pat_in(:)
|
|
|
|
contains
|
|
|
|
!-------------------------------------------------------------------
|
|
! compute vector product n1 of vectors v1 x v2
|
|
subroutine xprod(n1,v1,v2)
|
|
implicit none
|
|
|
|
real*8 n1(3), v1(3), v2(3)
|
|
|
|
n1(1) = v1(2)*v2(3) - v1(3)*v2(2)
|
|
n1(2) = v1(3)*v2(1) - v1(1)*v2(3)
|
|
n1(3) = v1(1)*v2(2) - v1(2)*v2(1)
|
|
|
|
end subroutine
|
|
|
|
|
|
!-------------------------------------------------------------------
|
|
! compute scalar product of vectors v1 and v2:
|
|
real*8 function scalar(v1,v2,n)
|
|
implicit none
|
|
integer i, n
|
|
real*8 v1(*), v2(*)
|
|
|
|
scalar=0.d0
|
|
do i=1,n
|
|
scalar=scalar+v1(i)*v2(i)
|
|
enddo
|
|
|
|
end function
|
|
|
|
!-------------------------------------------------------------------
|
|
! compute norm of vector:
|
|
real*8 function norm(x,n)
|
|
implicit none
|
|
integer i, n
|
|
real*8 x(*)
|
|
|
|
norm=0.d0
|
|
do i=1,n
|
|
norm=norm+x(i)**2
|
|
enddo
|
|
norm=sqrt(norm)
|
|
|
|
end function
|
|
|
|
!-------------------------------------------------------------------
|
|
|
|
subroutine rot_trisec(tr,v1,v2,v3)
|
|
implicit none
|
|
|
|
real*8 tr(3),v1(3),v2(3),v3(3)
|
|
|
|
real*8 vrot(3)
|
|
real*8 rot_ax(3)
|
|
real*8 cos_phi,sin_phi
|
|
|
|
! evaluate cos(-phi) and sin(-phi), where phi is the angle between
|
|
! tr and (1,0,0)
|
|
cos_phi=tr(1)
|
|
sin_phi=dsqrt(tr(2)**2+tr(3)**2)
|
|
|
|
if (sin_phi.lt.1.0d-12) then
|
|
return
|
|
endif
|
|
|
|
! determine rotational axis
|
|
rot_ax(1) = 0.0d0
|
|
rot_ax(2) = tr(3)
|
|
rot_ax(3) = -tr(2)
|
|
! normalize
|
|
rot_ax=rot_ax/sin_phi
|
|
|
|
! now the rotation can be done using Rodrigues' rotation formula
|
|
! v'=v*cos(p) + (k x v)sin(p) + k (k*v) (1-cos(p))
|
|
! for v=tr k*v vanishes by construction:
|
|
|
|
! check that the rotation does what it should
|
|
call rodrigues(vrot,tr,rot_ax,cos_phi,sin_phi)
|
|
if (dsqrt(vrot(2)**2+vrot(3)**2).gt.1.0d-12) then
|
|
write(6,*) "ERROR: BROKEN TRISECTOR"
|
|
stop
|
|
endif
|
|
tr=vrot
|
|
|
|
call rodrigues(vrot,v1,rot_ax,cos_phi,sin_phi)
|
|
v1=vrot
|
|
call rodrigues(vrot,v2,rot_ax,cos_phi,sin_phi)
|
|
v2=vrot
|
|
call rodrigues(vrot,v3,rot_ax,cos_phi,sin_phi)
|
|
v3=vrot
|
|
|
|
|
|
end subroutine
|
|
|
|
!-------------------------------------------------------------------
|
|
|
|
subroutine rodrigues(vrot,v,axis,cos_phi,sin_phi)
|
|
implicit none
|
|
real*8 vrot(3),v(3),axis(3)
|
|
real*8 cos_phi,sin_phi
|
|
|
|
real*8 ortho(3)
|
|
|
|
call xprod(ortho,axis,v)
|
|
vrot = v*cos_phi + ortho*sin_phi+axis*scalar(axis,v,3)*(1-cos_phi)
|
|
|
|
end subroutine
|
|
|
|
!-------------------------------------------------------------------
|
|
|
|
real*8 function get_ang(v,xaxis,yaxis)
|
|
implicit none
|
|
! get normalized [0:2pi) angle from vectors in the yz plane
|
|
real*8 v(3),xaxis(3),yaxis(3)
|
|
|
|
real*8 phi
|
|
|
|
real*8 pi
|
|
parameter (pi=3.141592653589793d0)
|
|
|
|
phi=atan2(scalar(yaxis,v,3),scalar(xaxis,v,3))
|
|
if (phi.lt.0.0d0) then
|
|
phi=2*pi+phi
|
|
endif
|
|
get_ang=phi
|
|
|
|
end function
|
|
|
|
end subroutine cart2int
|
|
subroutine genANN_ctrans(pat_in)
|
|
implicit none
|
|
|
|
include 'nnparams.incl'
|
|
include 'JTmod.incl'
|
|
|
|
real*8 pat_in(maxnin)
|
|
|
|
real*8 raw_in(maxnin),off_in(maxnin),ptrans_in(7)
|
|
real*8 r0
|
|
real*8 a,b,xs,ys,xb,yb
|
|
|
|
integer k
|
|
|
|
off_in(1:7)=pat_in(1:7)
|
|
r0=offsets(1)
|
|
|
|
! transform primitives
|
|
! recover raw distances from offset coords
|
|
do k=1,3
|
|
raw_in(k)=off_in(k)+offsets(1)
|
|
enddo
|
|
|
|
do k=1,3
|
|
ptrans_in(k)=off_in(k)
|
|
enddo
|
|
|
|
! rescale ONO angles
|
|
ptrans_in(4)=deg2rad*off_in(4)
|
|
ptrans_in(5)=deg2rad*off_in(5)
|
|
ptrans_in(6)=deg2rad*off_in(6)
|
|
! rescale umbrella
|
|
ptrans_in(7)=off_in(7)*deg2rad
|
|
|
|
! compute symmetry coordinates
|
|
|
|
! A (breathing)
|
|
a=(ptrans_in(1)+ptrans_in(2)+ptrans_in(3))/dsqrt(3.0d0)
|
|
! ES
|
|
call prim2emode(ptrans_in(1:3),xs,ys)
|
|
! EB
|
|
call prim2emode(ptrans_in(4:6),xb,yb)
|
|
! B (umbrella)
|
|
b=ptrans_in(7)
|
|
|
|
! overwrite input with output
|
|
|
|
pat_in(pat_index(1))=a ! 1
|
|
pat_in(pat_index(2))=xs
|
|
pat_in(pat_index(3))=ys
|
|
pat_in(pat_index(4))=xb
|
|
pat_in(pat_index(5))=yb
|
|
pat_in(pat_index(6))=b
|
|
! totally symmetric monomials
|
|
pat_in(pat_index(7))=xs**2 + ys**2 ! 2
|
|
pat_in(pat_index(8))=xb**2 + yb**2 ! 3
|
|
pat_in(pat_index(9))=b**2 ! 9
|
|
pat_in(pat_index(10))=xs*xb+ys*yb ! 4
|
|
! S^3, B^3
|
|
pat_in(pat_index(11))=xs*(xs**2-3*ys**2) ! 5
|
|
pat_in(pat_index(12))=xb*(xb**2-3*yb**2) ! 6
|
|
! S^2 B, S B^2
|
|
pat_in(pat_index(13))=xb*(xs**2-ys**2) - 2*yb*xs*ys ! 7
|
|
pat_in(pat_index(14))=xs*(xb**2-yb**2) - 2*ys*xb*yb ! 8
|
|
|
|
do k=11,14
|
|
pat_in(pat_index(k))=tanh(0.1d0*pat_in(pat_index(k)))*10.0d0
|
|
enddo
|
|
|
|
|
|
|
|
|
|
end subroutine
|
|
subroutine prim2emode(prim,ex,ey)
|
|
implicit none
|
|
! Takes a 2D-vector prim and returns the degenerate modes x and y
|
|
! following our standard conventions.
|
|
|
|
real*8 prim(3),ex,ey
|
|
|
|
ex=(2.0d0*prim(1)-prim(2)-prim(3))/dsqrt(6.0d0)
|
|
ey=(prim(2)-prim(3))/dsqrt(2.0d0)
|
|
|
|
end
|
|
end module ctrans_mod
|
|
|