-
Notifications
You must be signed in to change notification settings - Fork 17
/
ddcosmo.f90
1808 lines (1805 loc) · 53.4 KB
/
ddcosmo.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
module ddcosmo
implicit none
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
!
! 888 888 .d8888b. .d88888b. .d8888b. 888b d888 .d88888b.
! 888 888 d88P Y88b d88P" "Y88b d88P Y88b 8888b d8888 d88P" "Y88b
! 888 888 888 888 888 888 Y88b. 88888b.d88888 888 888
! .d88888 .d88888 888 888 888 "Y888b. 888Y88888P888 888 888
! d88" 888 d88" 888 888 888 888 "Y88b. 888 Y888P 888 888 888
! 888 888 888 888 888 888 888 888 "888 888 Y8P 888 888 888
! Y88b 888 Y88b 888 Y88b d88P Y88b. .d88P Y88b d88P 888 " 888 Y88b. .d88P
! "Y88888 "Y88888 "Y8888P" "Y88888P" "Y8888P" 888 888 "Y88888P"
!
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! COPYRIGHT (C) 2015 by Filippo Lipparini, Benjamin Stamm, Paolo Gatto !
! Eric Cancès, Yvon Maday, Jean-Philip Piquemal, Louis Lagardère and !
! Benedetta Mennucci. !
! ALL RIGHT RESERVED. !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
!
! A modular implementation of COSMO using a domain decomposition linear scaling
! strategy.
!
! This code is governed by the LGPL license and abiding by the rules of
! distribution of free software.
! This program is distributed in the hope that it will be useful, but
! WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
! or FITNESS FOR A PARTICULAR PURPOSE.
! See the GNU Lesser General Public License for more details.
!
! Users of this code are asked to include the following references in their
! publications:
!
! [1] E. Cancès, Y. Maday, B. Stamm
! "Domain decomposition for implicit solvation models"
! J. Chem. Phys. 139, 054111 (2013)
!
! [2] F. Lipparini, B. Stamm, E. Cancès, Y. Maday, B. Mennucci
! "Fast Domain Decomposition Algorithm for Continuum Solvation Models:
! Energy and First Derivatives"
! J. Chem. Theory Comput. 9, 3637–3648 (2013)
!
! Also, include one of the three following reference depending on whether you
! use this code in conjunction with a QM [3], Semiempirical [4] or Classical [5]
! description of the solute:
!
! [3] F. Lipparini, G. Scalmani, L. Lagardère, B. Stamm, E. Cancès, Y. Maday,
! J.-P. Piquemal, M. J. Frisch, B. Mennucci
! "Quantum, classical, and hybrid QM/MM calculations in solution: General
! implementation of the ddCOSMO linear scaling strategy"
! J. Chem. Phys. 141, 184108 (2014)
! (for quantum mechanical models)
!
! [4] F. Lipparini, L. Lagardère, G. Scalmani, B. Stamm, E. Cancès, Y. Maday,
! J.-P. Piquemal, M. J. Frisch, B. Mennucci
! "Quantum Calculations in Solution for Large to Very Large Molecules:
! A New Linear Scaling QM/Continuum Approach"
! J. Phys. Chem. Lett. 5, 953-958 (2014)
! (for semiempirical models)
!
! [5] F. Lipparini, L. Lagardère, C. Raynaud, B. Stamm, E. Cancès, B. Mennucci
! M. Schnieders, P. Ren, Y. Maday, J.-P. Piquemal
! "Polarizable Molecular Dynamics in a Polarizable Continuum Solvent"
! J. Chem. Theory Comput. 11, 623-634 (2015)
! (for classical models, including polarizable force fields
!
! The users of this code should also include the appropriate reference to the
! COSMO model. This distribution includes the routines to generate lebedev
! grids by D. Laikov and C. van Wuellen, as publicly available on CCL. If the
! routines are used, the following reference should also be included:
!
! [6] V.I. Lebedev, and D.N. Laikov
! "A quadrature formula for the sphere of the 131st
! algebraic order of accuracy"
! Doklady Mathematics, Vol. 59, No. 3, 1999, pp. 477-481.
!
! Written by Filippo Lipparini, October 2015.
!
!-------------------------------------------------------------------------------
!
!
! - numerical constants
!
integer, parameter :: ndiis=25, iout=6, nngmax=100
real*8, parameter :: zero=0.d0, pt5=0.5d0, one=1.d0, two=2.d0, four=4.d0
real*8, parameter :: se = -1.0d0
!
! - numerical constants explicitly computed
!
real*8 :: pi, sq2
!
! - quantities contained in control file
!
! iprint - printing flag
! nproc - number of openMP threads ; 0) no parallelization
! lmax - max angular momentum of spherical harmonics basis
! ngrid - desired number of Lebedev integration points
! iconv - threshold for iterative solver ( 10^-iconv )
! igrad - 1) compute forces ; 0) do not compute forces
! eps - dielectric constant of the solvent
! eta - regularization parameters
!
integer :: iprint, nproc, lmax, ngrid, iconv, igrad
real*8 :: eps, eta
!
! - other quantities
!
! nsph - number of spheres/atoms
! ncav - number of integration points on cavity's boundary
! nylm - number of basis functions, i.e., spherical harmonics
!
integer :: nsph, ncav, nylm
!
! - workspaces
!
integer, allocatable :: inl(:), nl(:)
real*8, allocatable :: rsph(:), csph(:,:), ccav(:,:)
real*8, allocatable :: w(:), grid(:,:), basis(:,:)
real*8, allocatable :: fact(:), facl(:), facs(:)
real*8, allocatable :: fi(:,:), ui(:,:), zi(:,:,:)
!
! - miscellanea
!
logical :: grad
!
! - subroutines & functions
! =========================
!
! ddinit - initialize data structure
! memfree - free data structure
! sprod - scalar product
! fsw - switch function
! dfsw - derivative of switch function
! ptcart - print routine
! prtsph - print routine
! intrhs - integrate spherical harmonics expansion
! ylmbas - compute spherical harmonics Y_l^m
! dbasis - compute derivatives of spherical harmonics
! polleg - compute Legendre polynomials
! trgev - service routine for computation of spherical harmonics
! intmlp - compute quantity needed for COSMO action
! wghpot - weigh potential at cavity points
! hsnorm - compute H-norm on unit sphere
! adjrhs - auxiliary routine for COSMO adjoint action
! header - print header
! fdoka - compute the first part of <S,L^(x)X>
! fdokb - compute the the second part of <S,L^(x)X>
! fdoga - compute the U^(x)\Phi contribution to <S,g^(x)>
! calcv - auxiliary routine for COSMO action
! ddmkxi - compute the xi intermediate
!
contains
!
!
!--------------------------------------------------------------------------------------------------
subroutine ddinit(n,x,y,z,rvdw)
implicit none
!
! allocate the various arrays needed for ddcosmo,
! assemble the cavity and the various associated geometrical quantities.
!
integer, intent(in) :: n
real*8, dimension(n), intent(in) :: x, y, z, rvdw
!
integer :: isph, jsph, i, ii, lnl, l, ind, m, igrid, inear, jnear
integer :: istatus
real*8 :: fac, fl, ffl, fnorm, d2, r2, v(3), vv, t, xt, swthr
!
real*8, allocatable :: vcos(:), vsin(:), vplm(:)
integer, parameter :: nllg=32
!
integer, dimension(nllg) :: ng0
data ng0/6,14,26,38,50,74,86,110,146,170,194,230,266,302,350,434,590,770,974, &
1202,1454,1730,2030,2354,2702,3074,3470,3890,4334,4802,5294,5810/
!
! openMP parallelization:
!
if (nproc.eq.0) nproc = 1
!$ call omp_set_num_threads(nproc)
pi = four*atan(one)
sq2 = sqrt(two)
nsph = n
!
! choose the lebedev grid with number of points closest to ngrid:
!
igrid = 0
inear = 100000
do i = 1, nllg
jnear = iabs(ng0(i)-ngrid)
if (jnear.lt.inear) then
inear = jnear
igrid = i
end if
end do
ngrid = ng0(igrid)
!
! print a nice header:
!
call header
!
! allocate:
!
grad = igrad.ne.0
nylm = (lmax+1)*(lmax+1)
allocate( rsph(nsph), csph(3,nsph), w(ngrid), grid(3,ngrid), basis(nylm,ngrid), &
inl(nsph+1), nl(nsph*nngmax), fi(ngrid,nsph), ui(ngrid,nsph), &
fact(max(2*lmax+1,2)), facl(nylm), facs(nylm) , stat=istatus )
if ( istatus.ne.0 ) then
write(*,*)'ddinit : [1] allocation failed !'
stop
endif
if (grad) allocate( zi(3,ngrid,nsph), stat=istatus)
if ( istatus.ne.0 ) then
write(*,*)'ddinit : [1] allocation failed !'
stop
endif
!
! precompute various quantities (diagonal blocks of ddcosmo matrix,
! normalization factors for spherical harmonics, factorials)
!
fact(1) = one
fact(2) = one
do i = 3, 2*lmax + 1
fact(i) = dble(i-1)*fact(i-1)
end do
!
do l = 0, lmax
ind = l*l + l + 1
fl = (two*dble(l) + one)/(four*pi)
ffl = sqrt(fl)
facl(ind-l:ind+l) = fl
facs(ind) = ffl
do m = 1, l
fnorm = sq2*ffl*sqrt(fact(l-m+1)/fact(l+m+1))
if (mod(m,2).eq.1) fnorm = -fnorm
facs(ind+m) = fnorm
facs(ind-m) = fnorm
end do
end do
!
! set the centers and radii of the spheres:
!
csph(1,:) = x
csph(2,:) = y
csph(3,:) = z
rsph = rvdw
!
! load a lebedev grid:
!
call llgrid(ngrid,w,grid)
!
! build a basis of spherical harmonics at the gridpoints:
!
allocate (vplm(nylm),vcos(lmax+1),vsin(lmax+1), stat=istatus)
if ( istatus.ne.0 ) then
write(*,*)'ddinit : [2] allocation failed !'
stop
endif
!$omp parallel do default(shared) private(i,vplm,vcos,vsin)
do i = 1, ngrid
call ylmbas(grid(:,i),basis(:,i),vplm,vcos,vsin)
end do
!$omp end parallel do
deallocate (vplm,vcos,vsin, stat=istatus)
if ( istatus.ne.0 ) then
write(*,*)'ddinit : [1] deallocation failed!'
stop
endif
!
if (iprint.ge.4) then
call prtsph('facs',1,0,facs)
call prtsph('facl',1,0,facs)
call prtsph('basis',ngrid,0,basis)
call ptcart('grid',3,0,grid)
call ptcart('weights',1,0,w)
end if
!
! build neighbors list (CSR format)
! =================================
!
! \\ jsph |
! isph \\ | 1 2 3 4 5 6
! -----------------------------------
! 1 | x x x
! 2 | x x x x
! 3 | x x x
! 4 | x x x x
! 5 | x x x
! 6 | x x x
!
!
! inl = 1, 4, 8, 11, 15, 18,21 pointer to 1st neighbor
!
! | | | | | |
! v V V V V V
!
! 1| 2| 3| 4| 5| 6| 7| 8| 9|10|11|12|13|14|15|16|17|18|19|20
!
! nl = 2, 4, 5, 1, 3, 5, 6, 2, 4, 6, 1, 3, 5, 6, 1, 2, 4, 2, 3, 4 neighbors list
!
!
! index of nl
!
ii = 1
lnl = 0
do isph = 1, nsph
inl(isph) = lnl + 1
do jsph = 1, nsph
if (isph.ne.jsph) then
d2 = (csph(1,isph) - csph(1,jsph))**2 + (csph(2,isph) - csph(2,jsph))**2 + (csph(3,isph) - csph(3,jsph))**2
r2 = (rsph(isph) + rsph(jsph))**2
if (d2.le.r2) then
nl(ii) = jsph
ii = ii + 1
lnl = lnl + 1
end if
end if
end do
end do
inl(nsph+1) = lnl+1
!
1000 format(t3,'neighbours of sphere ',i6)
1010 format(t5,12i6)
!
if (iprint.ge.4) then
write(iout,*) ' inl:'
write(iout,'(10i6)') inl(1:nsph+1)
write(iout,*)
do isph = 1, nsph
write(iout,1000) isph
write(iout,1010) nl(inl(isph):inl(isph+1)-1)
end do
write(iout,*)
end if
!
!-----------------------------------------------------------------------
! Define :
!
! N_i = list of neighbors of i-sphere [ excluding i-sphere ]
!
! | r_i + rho_i s_n - r_j |
! t_n^ij = -------------------------
! rho_j
!
! fi(n,i) = sum \chi( t_n^ij )
! j \in N_i
!
! Notice that the derivative of fi(n,i) wrt to r_k is (possibly) nonzero
! when either k = i, or k \in N_j .
!
!-----------------------------------------------------------------------
!
! build arrays fi, ui, zi
! =======================
fi = zero
ui = zero
if (grad) zi = zero
!
!$omp parallel do default(shared) private(isph,i,ii,jsph,v,vv,t,xt,swthr,fac)
!
! loop over spheres
do isph = 1,nsph
!
! loop over integration points
do i = 1,ngrid
!
! loop over neighbors of i-sphere
do ii = inl(isph),inl(isph+1)-1
!
! neighbor's number
jsph = nl(ii)
!
! compute t_n^ij
v(:) = csph(:,isph) + rsph(isph)*grid(:,i) - csph(:,jsph)
vv = sqrt(dot_product(v,v))
t = vv/rsph(jsph)
!
! compute \chi( t_n^ij )
xt = fsw( t, se, eta )
!
! upper bound of switch region
swthr = one + (se + 1.d0)*eta / 2.d0
!
! t_n^ij belongs to switch region
if ( grad .and. ( t.lt.swthr .and. t.gt.swthr-eta ) ) then
!
fac = dfsw( t, se, eta ) / rsph(jsph)
!
! accumulate for zi
! -----------------
zi(:,i,isph) = zi(:,i,isph) + fac*v(:)/vv
!
endif
!
! accumulate for fi
! -----------------
fi(i,isph) = fi(i,isph) + xt
!
enddo
!
! compute ui
! ----------
if ( fi(i,isph).le.one ) ui(i,isph) = one - fi(i,isph)
!
enddo
enddo
!$omp end parallel do
!
if (iprint.ge.4) then
call ptcart('fi',nsph,0,fi)
call ptcart('ui',nsph,0,ui)
end if
!
! build cavity array
! ==================
!
! initialize number of cavity points
ncav=0
!
! loop over spheres
do isph = 1,nsph
!
! loop over integration points
do i = 1,ngrid
!
! positive contribution from integration point
if ( ui(i,isph).gt.zero ) then
!
! accumulate
ncav = ncav + 1
!
endif
enddo
enddo
!
! allocate cavity array
allocate( ccav(3,ncav) , stat=istatus )
if ( istatus .ne. 0 ) then
write(*,*)'ddinit : [3] allocation failed!'
stop
endif
!
!
! initialize cavity array index
ii = 0
!
! loop over spheres
do isph = 1,nsph
!
! loop over integration points
do i = 1,ngrid
!
! positive contribution from integration point
if ( ui(i,isph).gt.zero ) then
!
! advance cavity array index
ii = ii + 1
!
! store point
ccav(:,ii) = csph(:,isph) + rsph(isph)*grid(:,i)
!
endif
enddo
enddo
!
1100 format(t3,i8,3f14.6)
!
if (iprint.ge.4) then
write(iout,*) ' external cavity points:'
do ii = 1, ncav
write(iout,1100) ii, ccav(:,ii)
end do
write(iout,*)
end if
!
return
!
!
end subroutine ddinit
!---------------------------------------------------------------------------------
!
!
!
!
!
!---------------------------------------------------------------------------------
subroutine memfree
!
implicit none
integer :: istatus, istatus0
!
!---------------------------------------------------------------------------------
!
! initialize deallocation flags
istatus0 = 0
istatus = 0
!
! deallocate the arrays
if( allocated(rsph) ) deallocate( rsph , stat=istatus ) ; istatus0 = istatus0 + istatus
if( allocated(csph) ) deallocate( csph , stat=istatus ) ; istatus0 = istatus0 + istatus
if( allocated(ccav) ) deallocate( ccav , stat=istatus ) ; istatus0 = istatus0 + istatus
if( allocated(w) ) deallocate( w , stat=istatus ) ; istatus0 = istatus0 + istatus
if( allocated(grid) ) deallocate( grid , stat=istatus ) ; istatus0 = istatus0 + istatus
if( allocated(basis) ) deallocate( basis , stat=istatus ) ; istatus0 = istatus0 + istatus
if( allocated(inl) ) deallocate( inl , stat=istatus ) ; istatus0 = istatus0 + istatus
if( allocated(nl) ) deallocate( nl , stat=istatus ) ; istatus0 = istatus0 + istatus
if( allocated(fact) ) deallocate( fact , stat=istatus ) ; istatus0 = istatus0 + istatus
if( allocated(facl) ) deallocate( facl , stat=istatus ) ; istatus0 = istatus0 + istatus
if( allocated(facs) ) deallocate( facs , stat=istatus ) ; istatus0 = istatus0 + istatus
if( allocated(ui) ) deallocate( ui , stat=istatus ) ; istatus0 = istatus0 + istatus
if( allocated(fi) ) deallocate( fi , stat=istatus ) ; istatus0 = istatus0 + istatus
if( allocated(zi) ) deallocate( zi , stat=istatus ) ; istatus0 = istatus0 + istatus
!
if ( istatus0.ne.0 ) then
write(*,*)'memfree : dallocation failed!'
stop
endif
!
end subroutine memfree
!---------------------------------------------------------------------------------
!
!
!
!
!
!---------------------------------------------------------------------------------
real*8 function sprod(n,u,v)
!
implicit none
integer, intent(in) :: n
real*8, dimension(n), intent(in) :: u, v
!
integer :: i
real*8 :: ss
!---------------------------------------------------------------------------------
!
! initialize
ss = zero
!
do i = 1, n
!
! accumulate
ss = ss + u(i)*v(i)
!
enddo
!
! redirect
sprod = ss
!
return
!
!
end function sprod
!------------------------------------------------------------------------------------------------
!
!
!
!
!------------------------------------------------------------------------------------------------
! Purpose : switch function (5th degree polynomial)
!------------------------------------------------------------------------------------------------
real*8 function fsw( t, s, eta )
!
implicit none
real*8, intent(in) :: t
real*8, intent(in) :: s
real*8, intent(in) :: eta
!
real*8 :: a, b, flow, x
real*8, parameter :: f6=6.0d0, f10=10.d0, f12=12.d0, f15=15.d0
!
!------------------------------------------------------------------------------------------------
!
! shift :
! s = 0 => t - eta/2 [ CENTERED ]
! s = 1 => t - eta [ EXTERIOR ]
! s = -1 => t [ INTERIOR ]
!
! apply shift
x = t - (s + 1.d0)*eta / 2.d0
!
! lower bound of switch region
flow = one - eta
!
! define switch function \chi
if ( x.ge.one ) then
!
fsw = zero
!
elseif ( x.le.flow ) then
!
fsw = one
!
else
!
a = f15*eta - f12
b = f10*eta*eta - f15*eta + f6
fsw = ( (x-one)*(x-one)*(one-x)*(f6*x*x + a*x + b) ) / (eta**5)
!
endif
!
!
end function fsw
!------------------------------------------------------------------------------------------------
!
!
!
!
!------------------------------------------------------------------------------------------------
! Purpose : first derivative of switch function
!------------------------------------------------------------------------------------------------
real*8 function dfsw( t, s, eta )
!
implicit none
real*8, intent(in) :: t
real*8, intent(in) :: s
real*8, intent(in) :: eta
!
real*8 flow, x
real*8, parameter :: f30=30.0d0
!
!------------------------------------------------------------------------------------------------
!
! shift :
! s = 0 => t - eta/2 [ CENTERED ]
! s = 1 => t - eta [ EXTERIOR ]
! s = -1 => t [ INTERIOR ]
!
! apply shift
x = t - (s + 1.d0)*eta / 2.d0
!
! lower bound of switch region
flow = one - eta
!
! define derivative of switch function \chi
if ( x.ge.one ) then
!
dfsw = zero
!
elseif ( x.le.flow ) then
!
dfsw = zero
!
else
!
dfsw = f30*(one-x)*(x-one)*(x-one+eta)*(x-one+eta)/(eta**5)
!
endif
!
!
end function dfsw
!------------------------------------------------------------------------------------------------
!
!
!
!
!
!
!
!------------------------------------------------------------------------------------------------
! Purpose : dump an array (ngrid,ncol) or just a column.
!------------------------------------------------------------------------------------------------
subroutine ptcart( label, ncol, icol, x )
!
implicit none
!
character (len=*), intent(in) :: label
integer, intent(in) :: ncol, icol
real*8, dimension(ngrid,ncol), intent(in) :: x
!
integer :: ig, noff, nprt, ic, j
!
!------------------------------------------------------------------------------------------------
!
! print header :
if (ncol.eq.1) then
write (iout,'(3x,a,1x,"(column ",i4")")') label, icol
else
write (iout,'(3x,a)') label
endif
!
! print entries :
if (ncol.eq.1) then
do ig = 1, ngrid
write(iout,1000) ig, x(ig,1)
enddo
!
else
noff = mod (ncol,5)
nprt = max(ncol - noff,0)
do ic = 1, nprt, 5
write(iout,1010) (j, j = ic, ic+4)
do ig = 1, ngrid
write(iout,1020) ig, x(ig,ic:ic+4)
end do
end do
write (iout,1010) (j, j = nprt+1, nprt+noff)
do ig = 1, ngrid
write(iout,1020) ig, x(ig,nprt+1:nprt+noff)
end do
end if
!
1000 format(1x,i5,f14.8)
1010 format(6x,5i14)
1020 format(1x,i5,5f14.8)
!
return
!
!
end subroutine ptcart
!------------------------------------------------------------------------------------------------
!
!
!
!
!
!------------------------------------------------------------------------------------------------
! Purpose : dump an array (nylm,ncol) or just a column.
!------------------------------------------------------------------------------------------------
subroutine prtsph(label,ncol,icol,x)
!
implicit none
!
character (len=*), intent(in) :: label
integer, intent(in) :: ncol, icol
real*8, dimension(nylm,ncol), intent(in) :: x
!
integer :: l, m, ind, noff, nprt, ic, j
!
!------------------------------------------------------------------------------------------------
!
! print header :
if (ncol.eq.1) then
write (iout,'(3x,a,1x,"(column ",i4")")') label, icol
else
write (iout,'(3x,a)') label
endif
!
! print entries :
if (ncol.eq.1) then
do l = 0, lmax
ind = l*l + l + 1
do m = -l, l
write(iout,1000) l, m, x(ind+m,1)
end do
end do
!
else
noff = mod (ncol,5)
nprt = max(ncol - noff,0)
do ic = 1, nprt, 5
write(iout,1010) (j, j = ic, ic+4)
do l = 0, lmax
ind = l*l + l + 1
do m = -l, l
write(iout,1020) l, m, x(ind+m,ic:ic+4)
end do
end do
end do
write (iout,1010) (j, j = nprt+1, nprt+noff)
do l = 0, lmax
ind = l*l + l + 1
do m = -l, l
write(iout,1020) l, m, x(ind+m,nprt+1:nprt+noff)
end do
end do
end if
!
1000 format(1x,i3,i4,f14.8)
1010 format(8x,5i14)
1020 format(1x,i3,i4,5f14.8)
!
return
!
!
end subroutine prtsph
!------------------------------------------------------------------------------------------------
!
!
!
!
!------------------------------------------------------------------------------------------------
! Purpose : integrate against spherical harmonics
!------------------------------------------------------------------------------------------------
subroutine intrhs( isph, x, xlm )
!
implicit none
integer, intent(in) :: isph
real*8, dimension(ngrid), intent(in) :: x
real*8, dimension(nylm), intent(inout) :: xlm
!
integer :: ig
!
!------------------------------------------------------------------------------------------------
!
! initialize
xlm = zero
!
! accumulate over integration points
do ig = 1, ngrid
!
xlm = xlm + basis(:,ig)*w(ig)*x(ig)
!
enddo
!
! printing
if ( iprint.ge.5 ) then
!
call ptcart('pot',1,isph,x)
call prtsph('vlm',1,isph,xlm)
!
endif
!
return
!
!
end subroutine intrhs
!------------------------------------------------------------------------------------------------
!
!
!
!
!------------------------------------------------------------------------------------------------
! Purpose : compute spherical harmonics
!------------------------------------------------------------------------------------------------
subroutine ylmbas( x, basloc, vplm, vcos, vsin )
!
implicit none
real*8, dimension(3), intent(in) :: x
real*8, dimension(nylm), intent(out) :: basloc, vplm
real*8, dimension(lmax+1), intent(out) :: vcos, vsin
!
integer :: l, m, ind
real*8 :: cthe, sthe, cphi, sphi, plm
!
!------------------------------------------------------------------------------------------------
!
! get cos(\theta), sin(\theta), cos(\phi) and sin(\phi) from the cartesian
! coordinates of x.
!
! evaluate cos( theta ) ; sin( theta )
cthe = x(3)
sthe = sqrt(one - cthe*cthe)
!
! evalutate cos( phi ) ; sin( phi )
if (sthe.ne.zero) then
cphi = x(1)/sthe
sphi = x(2)/sthe
else
cphi = zero
sphi = zero
endif
!
! evaluate cos(m*phi) and sin(m*phi) arrays. notice that this is
! pointless if z = 1, as the only non vanishing terms will be the
! ones with m=0.
if(sthe.ne.zero) then
call trgev(cphi,sphi,vcos,vsin)
else
vcos = one
vsin = zero
endif
!
! evaluate the generalized legendre polynomials
call polleg(cthe,sthe,vplm)
!
! now build the spherical harmonics. we will distinguish m=0,
! m>0 and m<0:
do l = 0, lmax
ind = l**2 + l + 1
!
! m = 0
basloc(ind) = facs(ind)*vplm(ind)
!
do m = 1, l
!
plm = vplm(ind+m)
!
! m > 0
basloc(ind+m) = facs(ind+m)*plm*vcos(m+1)
!
! m < 0
basloc(ind-m) = facs(ind-m)*plm*vsin(m+1)
!
enddo
enddo
!
return
!
!
end subroutine ylmbas
!------------------------------------------------------------------------------------------------
!
!
!
!
!------------------------------------------------------------------------------------------------
! Purpose : compute first derivatives of spherical harmonics
!------------------------------------------------------------------------------------------------
subroutine dbasis( x, basloc, dbsloc, vplm, vcos, vsin )
!
implicit none
real*8, dimension(3), intent(in) :: x
real*8, dimension(nylm), intent(inout) :: basloc, vplm
real*8, dimension(3,nylm), intent(inout) :: dbsloc
real*8, dimension(lmax+1), intent(inout) :: vcos, vsin
!
integer :: l, m, ind, VC, VS
real*8 :: cthe, sthe, cphi, sphi, plm, fln, pp1, pm1, pp
real*8 :: et(3), ep(3)
!
!------------------------------------------------------------------------------------------------
!
! get cos(\theta), sin(\theta), cos(\phi) and sin(\phi) from the cartesian
! coordinates of x.
cthe = x(3)
sthe = sqrt(one - cthe*cthe)
!
! not ( NORTH or SOUTH pole )
if ( sthe.ne.zero ) then
cphi = x(1)/sthe
sphi = x(2)/sthe
!
! NORTH or SOUTH pole
else
cphi = one
sphi = zero
end if
!
! evaluate the derivatives of theta and phi:
!
et(1) = cthe*cphi
et(2) = cthe*sphi
et(3) = -sthe
! not ( NORTH or SOUTH pole )
if( sthe.ne.zero ) then
ep(1) = -sphi/sthe
ep(2) = cphi/sthe
ep(3) = zero
!
! NORTH or SOUTH pole
else
ep(1) = zero
ep(2) = one
ep(3) = zero
!
end if
!
! evaluate cos(m*phi) and sin(m*phi) arrays. notice that this is
! pointless if z = 1, as the only non vanishing terms will be the
! ones with m=0.
!
! not ( NORTH or SOUTH pole )
if ( sthe.ne.zero ) then
!
call trgev( cphi, sphi, vcos, vsin )
!
! NORTH or SOUTH pole
else
!
vcos = one
vsin = zero
!
end if
VC=zero
VS=cthe
!
! evaluate the generalized legendre polynomials.
!
call polleg( cthe, sthe, vplm )
!
! now build the spherical harmonics. we will distinguish m=0,
! m>0 and m<0:
!
basloc = zero
dbsloc = zero
do l = 0, lmax
ind = l*l + l + 1
! m = 0
fln = facs(ind)
basloc(ind) = fln*vplm(ind)
if (l.gt.0) then
dbsloc(:,ind) = fln*vplm(ind+1)*et(:)
else
dbsloc(:,ind) = zero
end if
!dir$ simd
do m = 1, l
fln = facs(ind+m)
plm = fln*vplm(ind+m)
pp1 = zero
if (m.lt.l) pp1 = -pt5*vplm(ind+m+1)
pm1 = pt5*(dble(l+m)*dble(l-m+1)*vplm(ind+m-1))
pp = pp1 + pm1
!