The Elk Code
afindtstep.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2020 Peter Elliott, J. K. Dewhurst and S. Sharma.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: afindtstep
8 ! !INTERFACE:
9 subroutine afindtstep
10 ! !USES:
11 use modmain
12 use modtddft
13 use modmpi
14 ! !DESCRIPTION:
15 ! Performs a time step of the macroscopic Maxwell equation and updates the
16 ! induced vector potential ${\bf A}(t)$. In practice, a more general damped
17 ! Proca equation is solved:
18 ! $$ a_0{\bf A}+a_1\dot{\bf A}+a_2\ddot{\bf A}=
19 ! \frac{4\pi c}{\Omega}{\bf J}, $$
20 ! where $\Omega$ is the unit cell volume, ${\bf J}$ is the total current
21 ! across the unit cell, and the parameters $a_i$, $i=0,1,2$ are stored in the
22 ! array {\tt afindpm}. This generalisation allows for both a mass and damping
23 ! term, however the default values of $a_0=a_1=0$ and $a_2=1$ recover the
24 ! physical Maxwell equation.
25 !
26 ! !REVISION HISTORY:
27 ! Created January 2020 (P. Elliott)
28 ! Added mass and damping terms, December 2022 (JKD)
29 !EOP
30 !BOC
31 implicit none
32 ! local variables
33 integer i
34 real(8) dt,t1,t2,t3
35 ! time step length
36 dt=times(itimes+1)-times(itimes)
37 ! add to the time derivative of the induced A-field
39 t2=dt/afindpm(2)
40 do i=1,3
41  t3=t1*jtot(i)-afindpm(1)*afindt(i,1)-afindpm(0)*afindt(i,0)
42  afindt(i,1)=afindt(i,1)+t2*t3
43 end do
44 ! add to the induced A-field
45 afindt(1:3,0)=afindt(1:3,0)+afindt(1:3,1)*dt
46 ! add to the total A-field
47 afieldt(1:3,itimes+1)=afieldt(1:3,itimes+1)+afindt(1:3,0)
48 ! write the induced A-field and its time derivative to file
49 if (mp_mpi) then
50  open(50,file='AFINDT.OUT',form='FORMATTED')
51  write(50,'(6G18.10)') afindt(:,:)
52  close(50)
53 end if
54 end subroutine
55 !EOC
56 
logical mp_mpi
Definition: modmpi.f90:17
real(8), dimension(3) jtot
Definition: modmain.f90:748
real(8) omega
Definition: modmain.f90:20
real(8), dimension(3, 0:1) afindt
Definition: modtddft.f90:67
real(8) solsc
Definition: modmain.f90:1253
Definition: modmpi.f90:6
real(8), dimension(:,:), allocatable afieldt
Definition: modtddft.f90:58
real(8), dimension(:), allocatable times
Definition: modtddft.f90:48
real(8), dimension(0:2) afindpm
Definition: modtddft.f90:69
integer itimes
Definition: modtddft.f90:46
real(8), parameter fourpi
Definition: modmain.f90:1234
subroutine afindtstep
Definition: afindtstep.f90:10