Line data Source code
1 : !--------------------------------------------------------------------------------
2 : ! Copyright (c) 2016 Peter Grünberg Institut, Forschungszentrum Jülich, Germany
3 : ! This file is part of FLEUR and available as free software under the conditions
4 : ! of the MIT license as expressed in the LICENSE file in more detail.
5 : !--------------------------------------------------------------------------------
6 :
7 : MODULE m_types_stars
8 : USE m_juDFT
9 : IMPLICIT NONE
10 : PRIVATE
11 : PUBLIC :: t_stars
12 : !The stars
13 : TYPE:: t_stars
14 : !max-length of star
15 : REAL :: gmax=0.0
16 : REAL :: gmaxz=0.0
17 : !no of 3d-stars
18 : INTEGER :: ng3
19 : !no of 2d-stars
20 : INTEGER :: ng2
21 : !dim of box
22 : INTEGER ::mx1
23 : INTEGER ::mx2
24 : INTEGER ::mx3
25 :
26 :
27 : INTEGER :: ng3_fft !number of stars in fft-box of size 2*rkmax
28 :
29 : !rep. g-vector of star
30 : INTEGER, ALLOCATABLE ::kv3(:, :)
31 : !length of star
32 : REAL, ALLOCATABLE ::sk3(:)
33 : !mapping of g-vectors to stars
34 : INTEGER, ALLOCATABLE ::ig(:, :, :)
35 : !No of g-vectors in star
36 : INTEGER, ALLOCATABLE ::nstr(:)
37 : !rep. g-vector of 2D-star
38 : INTEGER, ALLOCATABLE ::kv2(:, :)
39 : !length of 2D-star
40 : REAL, ALLOCATABLE ::sk2(:)
41 : !No of g-vecs in 2D-star
42 : INTEGER, ALLOCATABLE ::nstr2(:)
43 : !mapping of
44 : INTEGER, ALLOCATABLE ::i2g(:,:) !map g_x,g_y to 2D-star
45 : INTEGER, ALLOCATABLE ::ig2(:) !map 3D-star to 2D-star
46 : INTEGER, ALLOCATABLE :: igvac(:,:) ! map 2D-star and g_z to 3D-star
47 :
48 : !
49 : REAL, ALLOCATABLE:: phi2(:) !for oneD code, currently not in use
50 : !phase factor of g-vector
51 : COMPLEX, ALLOCATABLE ::rgphs(:, :, :)
52 : !phase factor of 2d-g vectors
53 : COMPLEX, ALLOCATABLE ::r2gphs(:, :)
54 :
55 : COMPLEX, ALLOCATABLE :: ustep(:)
56 : REAL, ALLOCATABLE :: ufft(:)
57 :
58 : REAL, ALLOCATABLE :: ab3(:)
59 :
60 : ! q-shifted stuff
61 : REAL :: center(3) = [0.0,0.0,0.0]
62 : REAL, ALLOCATABLE :: gq(:, :)
63 : REAL, ALLOCATABLE :: gq2(:, :)
64 : COMPLEX, ALLOCATABLE :: ufft1(:)
65 : CONTAINS
66 : PROCEDURE :: mpi_bc=>mpi_bc_stars
67 : PROCEDURE :: init=>init_stars
68 : PROCEDURE :: dim=>dim_stars
69 : PROCEDURE :: map_2nd_vac
70 : PROCEDURE :: reset_stars
71 : END TYPE t_stars
72 : CONTAINS
73 160 : SUBROUTINE mpi_bc_stars(this,mpi_comm,irank)
74 : USE m_mpi_bc_tool
75 : CLASS(t_stars),INTENT(INOUT)::this
76 : INTEGER,INTENT(IN):: mpi_comm
77 : INTEGER,INTENT(IN),OPTIONAL::irank
78 : INTEGER ::rank,myrank,ierr,n
79 160 : IF (PRESENT(irank)) THEN
80 0 : rank=irank
81 : ELSE
82 160 : rank=0
83 : END IF
84 :
85 160 : CALL mpi_bc(this%gmax,rank,mpi_comm)
86 160 : CALL mpi_bc(this%gmaxz,rank,mpi_comm)
87 160 : CALL mpi_bc(this%ng3,rank,mpi_comm)
88 160 : CALL mpi_bc(this%ng2,rank,mpi_comm)
89 160 : CALL mpi_bc(this%mx1,rank,mpi_comm)
90 160 : CALL mpi_bc(this%mx2,rank,mpi_comm)
91 160 : CALL mpi_bc(this%mx3,rank,mpi_comm)
92 : !CALL mpi_bc(this%kimax,rank,mpi_comm)
93 : !CALL mpi_bc(this%kimax2,rank,mpi_comm)
94 : !CALL mpi_bc(this%kq1_fft,rank,mpi_comm)
95 : !CALL mpi_bc(this%kq2_fft,rank,mpi_comm)
96 : !CALL mpi_bc(this%kq3_fft,rank,mpi_comm)
97 : !CALL mpi_bc(this%kmxq_fft ,rank,mpi_comm)
98 : !CALL mpi_bc(this%igq_fft,rank,mpi_comm)
99 : !CALL mpi_bc(this%igq2_fft,rank,mpi_comm)
100 : !CALL mpi_bc(this%kxc1_fft,rank,mpi_comm)
101 : !CALL mpi_bc(this%kxc2_fft,rank,mpi_comm)
102 : !CALL mpi_bc(this%kxc3_fft,rank,mpi_comm)
103 160 : CALL mpi_bc(this%ng3_fft,rank,mpi_comm)
104 : !CALL mpi_bc(this%kmxxc_fft,rank,mpi_comm)
105 : !CALL mpi_bc(this%nxc3_fft ,rank,mpi_comm)
106 160 : CALL mpi_bc(this%kv3,rank,mpi_comm)
107 160 : CALL mpi_bc(this%sk3,rank,mpi_comm)
108 160 : CALL mpi_bc(this%ab3,rank,mpi_comm)
109 160 : CALL mpi_bc(this%ig,rank,mpi_comm)
110 160 : CALL mpi_bc(this%nstr,rank,mpi_comm)
111 160 : CALL mpi_bc(this%kv2,rank,mpi_comm)
112 160 : CALL mpi_bc(this%sk2,rank,mpi_comm)
113 160 : CALL mpi_bc(this%nstr2,rank,mpi_comm)
114 160 : CALL mpi_bc(this%ig2,rank,mpi_comm)
115 160 : CALL mpi_bc(this%i2g,rank,mpi_comm)
116 160 : CALL mpi_bc(this%phi2,rank,mpi_comm)
117 160 : CALL mpi_bc(this%rgphs,rank,mpi_comm)
118 160 : CALL mpi_bc(this%r2gphs,rank,mpi_comm)
119 :
120 160 : CALL mpi_bc(this%ustep,rank,mpi_comm)
121 160 : CALL mpi_bc(this%ufft,rank,mpi_comm)
122 :
123 160 : CALL mpi_bc(this%center(1),rank,mpi_comm)
124 160 : CALL mpi_bc(this%center(2),rank,mpi_comm)
125 160 : CALL mpi_bc(this%center(3),rank,mpi_comm)
126 160 : CALL mpi_bc(this%gq,rank,mpi_comm)
127 160 : CALL mpi_bc(this%gq2,rank,mpi_comm)
128 160 : CALL mpi_bc(this%ufft1,rank,mpi_comm)
129 :
130 160 : END SUBROUTINE mpi_bc_stars
131 :
132 80 : subroutine init_stars(stars,cell,sym,film,rkmax,qvec)
133 : USE m_spgrot
134 : USE m_types_cell
135 : USE m_types_sym
136 : USE m_sort
137 : CLASS(t_stars),INTENT(INOUT) :: stars
138 : TYPE(t_cell),INTENT(IN) :: cell
139 : TYPE(t_sym),INTENT(IN) :: sym
140 : LOGICAL,INTENT(IN) :: film
141 : REAL,INTENT(IN) :: rkmax
142 :
143 : REAL, OPTIONAL, INTENT(IN) :: qvec(3)
144 :
145 : INTEGER :: k1,k2,k3,n,n1,k
146 : REAL :: s,g(3),gmax2,sq
147 80 : INTEGER :: kr(3,sym%nop),kv(3)
148 80 : COMPLEX :: phas(sym%nop)
149 : LOGICAL :: l_insph
150 : INTEGER,ALLOCATABLE :: index(:)
151 80 : REAL,ALLOCATABLE :: gsk3(:), sk3q(:), sk2q(:)
152 :
153 : !modify with gmaxz here!
154 80 : gmax2=stars%gmax**2
155 400 : allocate(gsk3(stars%ng3),index(stars%ng3))
156 :
157 400 : ALLOCATE(stars%rgphs(-stars%mx1:stars%mx1,-stars%mx2:stars%mx2,-stars%mx3:stars%mx3))
158 400 : ALLOCATE(stars%kv3(3,stars%ng3),stars%sk3(stars%ng3),stars%nstr(stars%ng3))
159 80 : IF (stars%gmaxz>0.0) ALLOCATE(stars%ab3(stars%ng3))
160 80 : IF (PRESENT(qvec)) ALLOCATE(stars%gq(3,stars%ng3),sk3q(stars%ng3))
161 400 : ALLOCATE(stars%ig(-stars%mx1:stars%mx1,-stars%mx2:stars%mx2,-stars%mx3:stars%mx3))
162 :
163 1673730 : stars%rgphs=0.0
164 1673730 : stars%ig=0
165 276308 : stars%nstr=0
166 :
167 80 : k=0
168 : !Generate 3D stars
169 2062 : x_dim: DO k1 = stars%mx1,-stars%mx1,-1
170 1982 : kv(1) = k1
171 53986 : y_dim: DO k2 = stars%mx2,-stars%mx2,-1
172 51924 : kv(2) = k2
173 1664882 : z_dim: DO k3 = stars%mx3,-stars%mx3,-1
174 1610976 : IF ( stars%ig(k1,k2,k3) .NE. 0 ) CYCLE z_dim ! belongs to another star
175 1276002 : kv(3) = k3
176 :
177 20416032 : g=matmul(kv,cell%bmat)
178 : !IF (PRESENT(qvec)) g = g + matmul(qvec,cell%bmat) !!!!!!!!!!
179 5104008 : s=dot_product(g,g)
180 : !if (s>gmax2) cycle z_dim !not in sphere !!!!!!!!!!!!!
181 1276002 : sq=s
182 1276002 : IF (PRESENT(qvec)) THEN
183 0 : g = g + matmul(qvec,cell%bmat)
184 0 : sq=dot_product(g,g)
185 : END IF
186 1276002 : l_insph = sq>gmax2
187 1276002 : IF (ABS(stars%gmaxz).GE.1e-8) l_insph = .NOT.((dot_product(g(:2),g(:2)).LE.gmax2).AND.(ABS(g(3)).LE.stars%gmaxz))
188 1276002 : if (l_insph) cycle z_dim !not in sphere
189 276228 : k=k+1
190 1104912 : stars%kv3(:,k)=kv
191 276228 : stars%sk3(k)=sqrt(s)
192 276228 : IF (stars%gmaxz>0.0) stars%ab3(k) = ABS(g(3))
193 :
194 276228 : IF (PRESENT(qvec)) THEN
195 0 : stars%gq(:,k)=g
196 0 : sk3q(k)=sqrt(sq)
197 0 : stars%center=qvec
198 : END IF
199 :
200 : ! secondary key for equal length stars
201 : gsk3(k) = (stars%mx1+stars%kv3(1,k)) +&
202 : & (stars%mx2+stars%kv3(2,k))*(2*stars%mx1+1) +&
203 276228 : & (stars%mx3+stars%kv3(3,k))*(2*stars%mx1+1)*(2*stars%mx2+1)
204 :
205 : !Now generate all equivalent g-vectors
206 : CALL spgrot(sym%nop,sym%symor,sym%mrot,sym%tau,sym%invtab,stars%kv3(:,k),&
207 276228 : kr)
208 1030612 : DO n = 1,sym%nop
209 978688 : stars%ig(kr(1,n),kr(2,n),kr(3,n))=k
210 : ENDDO
211 : ENDDO z_dim
212 : ENDDO y_dim
213 : ENDDO x_dim
214 80 : if (k.ne.stars%ng3) call judft_error("BUG inconsistency in star setup")
215 :
216 : !sort for increasing length sk3
217 80 : CALL sort(index,stars%sk3,gsk3)
218 2209984 : stars%kv3=stars%kv3(:,index)
219 552696 : stars%sk3=stars%sk3(index)
220 80 : IF (stars%gmaxz>0.0) stars%ab3=stars%ab3(index)
221 80 : IF (PRESENT(qvec)) stars%gq=stars%gq(:,index)
222 0 : IF (PRESENT(qvec)) stars%sk3=sk3q(index)
223 : ! set up the pointers and phases for 3d stars
224 276308 : DO k = 1,stars%ng3
225 : CALL spgrot(sym%nop,sym%symor,sym%mrot,sym%tau,sym%invtab,stars%kv3(:,k),&
226 276228 : kr,phas)
227 978768 : symloop: DO n = 1,sym%nop
228 702460 : stars%ig(kr(1,n),kr(2,n),kr(3,n))=k
229 702460 : stars%rgphs(kr(1,n),kr(2,n),kr(3,n))=stars%rgphs(kr(1,n),kr(2,n),kr(3,n))+phas(n)
230 3827479 : DO n1 = 1,n-1
231 4792799 : if (all(kr(:,n)==kr(:,n1))) cycle symloop !This g-vector we have already
232 : ENDDO
233 978688 : stars%nstr(k)=stars%nstr(k)+1
234 : ENDDO symloop
235 : ENDDO
236 : !Adjust phases
237 80 : if (sym%symor) THEN
238 1530538 : stars%rgphs=1.0
239 : ELSE
240 194 : DO k1 = stars%mx1,-stars%mx1,-1
241 5425 : DO k2 = stars%mx2,-stars%mx2,-1
242 143617 : DO k3 = stars%mx3,-stars%mx3,-1
243 138199 : IF ( stars%ig(k1,k2,k3)==0 ) CYCLE
244 143430 : stars%rgphs(k1,k2,k3)=stars%rgphs(k1,k2,k3)*stars%nstr(stars%ig(k1,k2,k3))/sym%nop
245 : enddo
246 : ENDDO
247 : ENDDO
248 : ENDIF
249 : !count number of stars in 2*rkmax (stars are ordered)
250 : associate(i=>stars%ng3_fft)
251 192066 : DO i=stars%ng3,1,-1
252 192066 : IF ( stars%sk3(i).LE.2.0*rkmax ) EXIT
253 : ENDDO
254 : end associate
255 :
256 287 : if (.not.film) return
257 : !
258 : ! Now the same for the 2D stars
259 : !
260 44 : ALLOCATE(stars%r2gphs(-stars%mx1:stars%mx1,-stars%mx2:stars%mx2))
261 77 : ALLOCATE(stars%kv2(2,stars%ng2),stars%sk2(stars%ng2),stars%nstr2(stars%ng2))
262 11 : IF (PRESENT(qvec)) ALLOCATE(stars%gq2(2,stars%ng3),sk2q(stars%ng2))
263 44 : ALLOCATE(stars%i2g(-stars%mx1:stars%mx1,-stars%mx2:stars%mx2))
264 33 : ALLOCATE(stars%ig2(stars%ng3))
265 44 : ALLOCATE(stars%igvac(stars%ng2,-stars%mx3:stars%mx3))
266 7987 : stars%r2gphs=0.0
267 7987 : stars%i2g=0
268 2035 : stars%nstr2=0
269 : kv(3)=0
270 11 : k=0
271 : !Generate 2D stars
272 296 : x_dim2: DO k1 = stars%mx1,-stars%mx1,-1
273 285 : kv(1) = k1
274 7983 : y_dim2: DO k2 = stars%mx2,-stars%mx2,-1
275 7687 : kv(2) = k2
276 7687 : IF ( stars%i2g(k1,k2) .NE. 0 ) CYCLE y_dim2 ! belongs to another star
277 51942 : g(:2)=matmul(kv(:2),cell%bmat(:2,:2))
278 14166 : s=dot_product(g(:2),g(:2))
279 4722 : sq=s
280 4722 : IF (PRESENT(qvec)) THEN
281 0 : g(:2) = g(:2) + matmul(qvec(:2),cell%bmat(:2,:2))
282 0 : sq=dot_product(g(:2),g(:2))
283 : END IF
284 4722 : if (sq>gmax2) cycle y_dim2 !not in sphere
285 2024 : k=k+1
286 6072 : stars%kv2(:,k)=kv(:2)
287 2024 : stars%sk2(k)=sqrt(s)
288 2024 : IF (PRESENT(qvec)) THEN
289 0 : stars%gq2(:2,k)=g(:2)
290 0 : sk2q(k)=sqrt(sq)
291 : END IF
292 : ! secondary key for equal length stars
293 2024 : gsk3(k) = (stars%mx1+stars%kv3(1,k)) +(stars%mx2+stars%kv3(2,k))*(2*stars%mx1+1)
294 : !Now generate all equivalent g-vectors
295 2024 : CALL spgrot(sym%nop2,sym%symor,sym%mrot(:2,:2,:),sym%tau(:2,:),sym%invtab,stars%kv2(:2,k),kr(:2,:))
296 7755 : DO n = 1,sym%nop2
297 7470 : stars%i2g(kr(1,n),kr(2,n))=k
298 : ENDDO
299 : ENDDO y_dim2
300 : ENDDO x_dim2
301 11 : if (k.ne.stars%ng2) call judft_error("BUG in init_stars: inconsistency in ng2")
302 : !sort for increasing length sk2
303 11 : CALL sort(index(:stars%ng2),stars%sk2,gsk3(:stars%ng2))
304 12155 : stars%kv2(:,:)=stars%kv2(:,index(:stars%ng2))
305 4081 : stars%sk2=stars%sk2(index(:stars%ng2))
306 11 : IF (PRESENT(qvec)) stars%gq2=stars%gq2(:,index(:stars%ng2))
307 0 : IF (PRESENT(qvec)) stars%sk2=sk2q(index(:stars%ng2))
308 : ! set up the pointers and phases for 2d stars
309 2035 : DO k = 1,stars%ng2
310 80074 : DO k3= stars%mx3,-stars%mx3,-1
311 80074 : stars%igvac(k,k3) = stars%ig(stars%kv2(1,k),stars%kv2(2,k),k3)
312 : ENDDO
313 2024 : CALL spgrot(sym%nop2,sym%symor,sym%mrot(:2,:2,:),sym%tau(:2,:),sym%invtab,stars%kv2(:2,k),kr(:2,:),phas)
314 7481 : symloop2: DO n = 1,sym%nop2
315 5446 : stars%i2g(kr(1,n),kr(2,n))=k
316 5446 : stars%r2gphs(kr(1,n),kr(2,n))=stars%r2gphs(kr(1,n),kr(2,n))+phas(n)
317 13652 : DO n1 = 1,n-1
318 15492 : if (all(kr(:2,n)==kr(:2,n1))) cycle symloop2 !This g-vector we have already
319 : ENDDO
320 7470 : stars%nstr2(k)=stars%nstr2(k)+1
321 : ENDDO symloop2
322 : ENDDO
323 29366 : DO k=1,stars%ng3
324 29366 : stars%ig2(k)=stars%i2g(stars%kv3(1,k),stars%kv3(2,k))
325 : ENDDO
326 : !Adjust phases
327 11 : IF (sym%symor) THEN
328 7987 : stars%r2gphs=1.0
329 : ELSE
330 0 : DO k1 = stars%mx1,-stars%mx1,-1
331 0 : DO k2 = stars%mx2,-stars%mx2,-1
332 0 : IF ( stars%i2g(k1,k2)==0 ) CYCLE
333 0 : stars%r2gphs(k1,k2)=stars%r2gphs(k1,k2)*stars%nstr2(stars%i2g(k1,k2))/sym%nop2
334 : ENDDO
335 : ENDDO
336 : ENDIF
337 287 : END SUBROUTINE init_stars
338 :
339 80 : subroutine dim_stars(stars,sym,cell,film,qvec)
340 : !! determine the key dimensions of the stars:
341 : !! mx1,mx2,mx3
342 : !! ng3, ng2
343 : USE m_spgrot
344 : USE m_types_cell
345 : USE m_types_sym
346 : USE m_boxdim
347 : CLASS(t_stars),INTENT(INOUT) :: stars
348 : TYPE(t_cell),INTENT(IN) :: cell
349 : TYPE(t_sym),INTENT(IN) :: sym
350 : LOGICAL,INTENT(IN) :: film
351 :
352 : REAL, OPTIONAL, INTENT(IN) :: qvec(3)
353 :
354 80 : INTEGER :: k1,k2,k3,n,kv(3),kr(3,sym%nop)
355 : REAL :: s,g(3),gmax2
356 : REAL :: arltv1,arltv2,arltv3
357 : LOGICAL :: l_insph
358 :
359 80 : gmax2=stars%gmax**2
360 80 : CALL boxdim(cell%bmat,arltv1,arltv2,arltv3)
361 :
362 80 : stars%mx1 = int(stars%gmax/arltv1) + 1
363 80 : stars%mx2 = int(stars%gmax/arltv2) + 1
364 80 : stars%mx3 = int(stars%gmax/arltv3) + 1
365 80 : IF (ABS(stars%gmaxz).GE.1e-8) stars%mx3 = int(stars%gmaxz/arltv3) + 1
366 400 : ALLOCATE(stars%ig(-stars%mx1:stars%mx1,-stars%mx2:stars%mx2,-stars%mx3:stars%mx3))
367 320 : ALLOCATE(stars%i2g(-stars%mx1:stars%mx1,-stars%mx2:stars%mx2))
368 :
369 53990 : stars%i2g=0
370 1673730 : stars%ig=0
371 80 : stars%ng2=0
372 80 : stars%ng3=0
373 2062 : x_dim: DO k1 = stars%mx1,-stars%mx1,-1
374 1982 : kv(1) = k1
375 53986 : y_dim: DO k2 = stars%mx2,-stars%mx2,-1
376 51924 : kv(2) = k2
377 51924 : kv(3) = 0
378 : !Check 2d-star
379 51924 : IF (stars%i2g(k1,k2)==0) THEN
380 614816 : g=matmul(kv,cell%bmat)
381 38426 : IF (PRESENT(qvec)) g = g + matmul(qvec,cell%bmat)
382 153704 : s=dot_product(g,g)
383 38426 : l_insph = .not.s>gmax2
384 38426 : IF (l_insph) THEN !in sphere
385 17758 : stars%ng2=stars%ng2+1
386 17758 : CALL spgrot(sym%nop2,sym%symor,sym%mrot,sym%tau,sym%invtab,kv,kr)
387 97574 : DO n = 1,sym%nop2
388 79816 : if (kr(3,n).ne.0) cycle
389 97574 : stars%i2g(kr(1,n),kr(2,n))=stars%ng2
390 : ENDDO
391 : ENDIF
392 : ENDIF
393 1664882 : z_dim: DO k3 = stars%mx3,-stars%mx3,-1
394 1610976 : IF ( stars%ig(k1,k2,k3) .NE. 0 ) CYCLE z_dim ! belongs to another star
395 1276002 : kv(3) = k3
396 20416032 : g=matmul(kv,cell%bmat)
397 1276002 : IF (PRESENT(qvec)) g = g + matmul(qvec,cell%bmat)
398 5104008 : s=dot_product(g,g)
399 1276002 : l_insph = .not.s>gmax2
400 1276002 : IF (ABS(stars%gmaxz).GE.1e-8) l_insph = (dot_product(g(:2),g(:2)).LE.gmax2).AND.(ABS(g(3)).LE.stars%gmaxz)
401 1276002 : if (.NOT.l_insph) cycle z_dim !not in sphere
402 276228 : stars%ng3=stars%ng3+1
403 276228 : CALL spgrot(sym%nop,sym%symor,sym%mrot,sym%tau,sym%invtab,kv(:),kr)
404 1030612 : DO n = 1,sym%nop
405 978688 : stars%ig(kr(1,n),kr(2,n),kr(3,n))=stars%ng3
406 : ENDDO
407 : ENDDO z_dim
408 : ENDDO y_dim
409 : ENDDO x_dim
410 :
411 80 : if (.not.film) stars%ng2=0
412 80 : DEALLOCATE(stars%ig)
413 80 : DEALLOCATE(stars%i2g)
414 :
415 :
416 80 : END SUBROUTINE dim_stars
417 :
418 126158 : subroutine map_2nd_vac(stars,vacuum,n2,n2_rot,phas)
419 : USE m_constants,ONLY: tpi_const
420 : USE m_types_vacuum
421 : CLASS(t_stars),INTENT(IN) :: stars
422 : TYPE(t_vacuum),INTENT(IN) :: vacuum
423 : INTEGER,INTENT(IN) :: n2
424 : INTEGER,INTENT(OUT) :: n2_rot
425 : COMPLEX,INTENT(OUT) :: phas
426 :
427 : INTEGER:: kr(2)
428 : REAL :: arg
429 :
430 1892370 : kr = matmul(stars%kv2(:,n2),vacuum%mrot2)
431 126158 : n2_rot = stars%i2g(kr(1),kr(2))
432 126158 : arg = tpi_const* ( stars%kv2(1,n2)*vacuum%tau2(1) + stars%kv2(2,n2)*vacuum%tau2(2) )
433 126158 : phas = cmplx(cos(arg),sin(arg))
434 :
435 126158 : end subroutine
436 :
437 0 : SUBROUTINE reset_stars(stars)
438 : CLASS(t_stars), INTENT(INOUT) :: stars
439 :
440 0 : stars%gmax=0.0
441 0 : stars%gmaxz=0.0
442 0 : stars%center = [0.0,0.0,0.0]
443 0 : IF (ALLOCATED(stars%kv3)) DEALLOCATE(stars%kv3)
444 0 : IF (ALLOCATED(stars%sk3)) DEALLOCATE(stars%sk3)
445 0 : IF (ALLOCATED(stars%ab3)) DEALLOCATE(stars%ab3)
446 0 : IF (ALLOCATED(stars%ig)) DEALLOCATE(stars%ig)
447 0 : IF (ALLOCATED(stars%nstr)) DEALLOCATE(stars%nstr)
448 0 : IF (ALLOCATED(stars%rgphs)) DEALLOCATE(stars%rgphs)
449 0 : IF (ALLOCATED(stars%ustep)) DEALLOCATE(stars%ustep)
450 0 : IF (ALLOCATED(stars%ufft)) DEALLOCATE(stars%ufft)
451 0 : IF (ALLOCATED(stars%gq)) DEALLOCATE(stars%gq)
452 0 : IF (ALLOCATED(stars%gq2)) DEALLOCATE(stars%gq2)
453 0 : IF (ALLOCATED(stars%ufft1)) DEALLOCATE(stars%ufft1)
454 0 : IF (ALLOCATED(stars%kv2)) DEALLOCATE(stars%kv2)
455 0 : IF (ALLOCATED(stars%sk2)) DEALLOCATE(stars%sk2)
456 0 : IF (ALLOCATED(stars%nstr2)) DEALLOCATE(stars%nstr2)
457 0 : IF (ALLOCATED(stars%i2g)) DEALLOCATE(stars%i2g)
458 0 : IF (ALLOCATED(stars%ig2)) DEALLOCATE(stars%ig2)
459 0 : IF (ALLOCATED(stars%igvac)) DEALLOCATE(stars%igvac)
460 0 : IF (ALLOCATED(stars%phi2)) DEALLOCATE(stars%phi2)
461 0 : IF (ALLOCATED(stars%r2gphs)) DEALLOCATE(stars%r2gphs)
462 0 : END SUBROUTINE
463 0 : END MODULE m_types_stars
|