![]() |
RFOF - RF Library for Orbit Following Codes
|
00001 00007 00008 module euitm_waves_interface 00009 00010 use RFOF_parameters 00011 use euitm_schemas 00012 00013 implicit none 00014 00015 contains 00016 00019 subroutine construct_dummy_euitm_waves(waves) 00020 00021 ! Output 00022 type(type_waves), intent(out) :: waves 00023 ! Local 00024 integer :: nfreq, nphi, nnphi, nion 00025 integer :: npsi, ntheta, max_npsi, max_ntheta 00026 integer :: jfreq, jnphi, jion, jpsi, jtheta 00027 00028 real(8) :: RFpower, EfieldNormalisation, freq, kperp 00029 00030 NAMELIST/input_wavefields/nfreq, nnphi, RFpower, & 00031 EfieldNormalisation, freq, nphi, kperp 00032 00033 open( io_channel_3872, FILE='input.rfof') 00034 read( io_channel_3872, input_wavefields) 00035 close(io_channel_3872) 00036 00037 nion = 1 00038 npsi = 1 00039 ntheta = 1 00040 00041 max_npsi = npsi 00042 max_ntheta = ntheta 00043 00044 allocate(waves%coherentwave(nfreq)) 00045 00046 !---- waves%global_param ----! 00047 do jfreq = 1, nfreq 00048 allocate(waves%coherentwave(jfreq)%global_param%ntor(nnphi)) 00049 allocate(waves%coherentwave(jfreq)%global_param%p_frac_ntor(nnphi)) 00050 allocate(waves%coherentwave(jfreq)%global_param%pow_i(nion)) 00051 allocate(waves%coherentwave(jfreq)%global_param%pow_ntor_i(nnphi,nion)) 00052 allocate(waves%coherentwave(jfreq)%global_param%pow_ntor_e(nnphi)) 00053 allocate(waves%coherentwave(jfreq)%global_param%cur_tor_ntor(nnphi)) 00054 00055 allocate(waves%coherentwave(jfreq)%global_param%name(1)) 00056 allocate(waves%coherentwave(jfreq)%global_param%type(1)) 00057 00058 waves%coherentwave(jfreq)%global_param%frequency = freq * (1d0+0.01d0*dble(jfreq-1)) 00059 waves%coherentwave(jfreq)%global_param%name = 'Antenna made up in RFOF' 00060 waves%coherentwave(jfreq)%global_param%type = 'IC' 00061 waves%coherentwave(jfreq)%global_param%ntor = nphi 00062 ! waves%coherentwave(jfreq)%global_param%f_assumption ! No use here! 00063 waves%coherentwave(jfreq)%global_param%power_tot = RFpower 00064 waves%coherentwave(jfreq)%global_param%pow_e = 0d0 00065 waves%coherentwave(jfreq)%global_param%cur_tor = 0d0 00066 ! waves%coherentwave(jfreq)%global_param%code_type ! No use here! 00067 ! waves%coherentwave(jfreq)%global_param%toroid_field... ! No use here! 00068 00069 do jion = 1, nion 00070 waves%coherentwave(jfreq)%global_param%pow_i(jion) = RFpower / real(nion) 00071 enddo 00072 00073 do jnphi = 1, nnphi 00074 waves%coherentwave(jfreq)%global_param%ntor( jnphi) = nphi * jnphi 00075 waves%coherentwave(jfreq)%global_param%pow_ntor_e( jnphi) = 0d0 00076 waves%coherentwave(jfreq)%global_param%p_frac_ntor( jnphi) = 1d0 / real(nnphi) 00077 waves%coherentwave(jfreq)%global_param%cur_tor_ntor( jnphi) = 0d0 00078 00079 do jion = 1, nion 00080 waves%coherentwave(jfreq)%global_param%pow_ntor_i(jnphi,jion) = RFpower / real(nnphi * nion) 00081 enddo 00082 enddo 00083 enddo 00084 !print *, "itm_waves: global_param: nnphi", waves%global_param%nntor 00085 !print *, "itm_waves: global_param: nphi", waves%global_param%ntor 00086 !print *, "itm_waves: global_param: ptot", waves%global_param%power_tot 00087 !print *, "itm_waves: global_param: j", waves%global_param%cur_tor_ntor 00088 !print *, "itm_waves: global_param: freq", waves%global_param%frequency 00089 00090 !---- waves%fullwave%local ----! 00091 do jfreq = 1, nfreq 00092 allocate( waves%coherentwave(jfreq)%fullwave%local%e_plus( nnphi, max_npsi,max_ntheta) ) 00093 allocate( waves%coherentwave(jfreq)%fullwave%local%e_minus( nnphi, max_npsi,max_ntheta) ) 00094 allocate( waves%coherentwave(jfreq)%fullwave%local%e_para( nnphi, max_npsi,max_ntheta) ) 00095 00096 allocate( waves%coherentwave(jfreq)%fullwave%local%e_plus_ph( nnphi, max_npsi,max_ntheta) ) 00097 allocate( waves%coherentwave(jfreq)%fullwave%local%e_minus_ph( nnphi, max_npsi,max_ntheta) ) 00098 allocate( waves%coherentwave(jfreq)%fullwave%local%e_para_ph( nnphi, max_npsi,max_ntheta) ) 00099 00100 do jnphi = 1, nnphi 00101 do jpsi = 1, max_npsi 00102 do jtheta = 1, max_ntheta 00103 waves%coherentwave(jfreq)%fullwave%local%e_plus( jnphi, jpsi,jtheta) = 0.5d0 00104 waves%coherentwave(jfreq)%fullwave%local%e_minus( jnphi, jpsi,jtheta) = 0.5d0 00105 waves%coherentwave(jfreq)%fullwave%local%e_para( jnphi, jpsi,jtheta) = 0d0 00106 waves%coherentwave(jfreq)%fullwave%local%e_plus_ph( jnphi, jpsi,jtheta) = 0d0 00107 waves%coherentwave(jfreq)%fullwave%local%e_minus_ph(jnphi, jpsi,jtheta) = 0d0 00108 waves%coherentwave(jfreq)%fullwave%local%e_para_ph( jnphi, jpsi,jtheta) = 0d0 00109 enddo 00110 enddo 00111 enddo 00112 enddo 00113 00114 !print *, "e+ = ",waves%fullwave%local%e_plus 00115 !print *, "e- = ",waves%fullwave%local%e_minus 00116 !print *, "ph e+ = ",waves%fullwave%local%e_plus_ph 00117 !print *, "ph e- = ",waves%fullwave%local%e_minus_ph 00118 00119 end subroutine construct_dummy_euitm_waves 00120 00121 !!$ !> Fill in a dummy waves CPO. 00122 !!$ !> Using version phase 4.08a, generated 21/04/2010 00123 !!$ subroutine construct_dummy_euitm_waves(waves) 00124 !!$ 00125 !!$ ! Output 00126 !!$ type(type_waves), intent(out) :: waves 00127 !!$ ! Local 00128 !!$ integer :: nfreq, nphi, nnphi, max_nnphi, nion 00129 !!$ integer :: npsi, ntheta, max_npsi, max_ntheta 00130 !!$ integer :: jfreq, jnphi, jion, jpsi, jtheta 00131 !!$ 00132 !!$ real(8) :: RFpower, EfieldNormalisation, freq, kperp 00133 !!$ 00134 !!$ NAMELIST/input_wavefields/nfreq, nnphi, RFpower, & 00135 !!$ EfieldNormalisation, freq, nphi, kperp 00136 !!$ 00137 !!$ open( io_channel_3872, FILE='input.rfof') 00138 !!$ read( io_channel_3872, input_wavefields) 00139 !!$ close(io_channel_3872) 00140 !!$ 00141 !!$ max_nnphi = nnphi 00142 !!$ 00143 !!$ nion = 1 00144 !!$ npsi = 1 00145 !!$ ntheta = 1 00146 !!$ 00147 !!$ max_npsi = npsi 00148 !!$ max_ntheta = ntheta 00149 !!$ 00150 !!$ !---- waves%global_param ----! 00151 !!$ 00152 !!$ !print *, "ok tryin itm waves stuff...allocate" 00153 !!$ 00154 !!$ allocate(waves%global_param%frequency(nfreq)) 00155 !!$ allocate(waves%global_param%name(nfreq)) 00156 !!$ allocate(waves%global_param%type(nfreq)) 00157 !!$ allocate(waves%global_param%nntor(nfreq)) 00158 !!$ allocate(waves%global_param%ntor(nfreq,max_nnphi)) 00159 !!$ ! ERROR in type description!!! allocate(waves%global_param%f_assumption(nion+1)) 00160 !!$ 00161 !!$ allocate(waves%global_param%power_tot(nfreq)) 00162 !!$ allocate(waves%global_param%p_frac_ntor(nfreq,max_nnphi)) 00163 !!$ allocate(waves%global_param%pow_i(nfreq,nion)) 00164 !!$ allocate(waves%global_param%pow_e(nfreq)) 00165 !!$ 00166 !!$ allocate(waves%global_param%pow_ntor_i(nfreq,nnphi,nion)) 00167 !!$ allocate(waves%global_param%pow_ntor_e(nfreq,nnphi)) 00168 !!$ allocate(waves%global_param%cur_tor(nfreq)) 00169 !!$ allocate(waves%global_param%cur_tor_ntor(nfreq,nnphi)) 00170 !!$ allocate(waves%global_param%code_type(nfreq)) 00171 !!$ ! Something beam tracing - who care :) !! allocate(waves%global_param%freq_point(nfreq)) 00172 !!$ 00173 !!$ do jfreq = 1, nfreq 00174 !!$ waves%global_param%frequency(jfreq) = freq * (1d0+0.01d0*dble(jfreq-1)) 00175 !!$ waves%global_param%name(jfreq) = 'No antenna name' 00176 !!$ waves%global_param%type(jfreq) = 'IC' 00177 !!$ waves%global_param%nntor(jfreq) = nnphi 00178 !!$ waves%global_param%power_tot(jfreq) = RFpower 00179 !!$ waves%global_param%pow_e(jfreq) = 0d0 00180 !!$ waves%global_param%cur_tor(jfreq) = 0d0 00181 !!$ waves%global_param%code_type(jfreq) = 0 00182 !!$ 00183 !!$ do jion = 1, nion 00184 !!$ waves%global_param%pow_i( jfreq,jion) = RFpower / real(nion) 00185 !!$ enddo 00186 !!$ 00187 !!$ do jnphi = 1, waves%global_param%nntor(jfreq) 00188 !!$ waves%global_param%ntor( jfreq,jnphi) = nphi * jnphi 00189 !!$ waves%global_param%pow_ntor_e( jfreq,jnphi) = 0d0 00190 !!$ waves%global_param%p_frac_ntor( jfreq,jnphi) = 1d0 / real(nnphi) 00191 !!$ waves%global_param%cur_tor_ntor( jfreq,jnphi) = 0d0 00192 !!$ 00193 !!$ do jion = 1, nion 00194 !!$ waves%global_param%pow_ntor_i(jfreq,jnphi,jion) = RFpower / real(nnphi * nion) 00195 !!$ enddo 00196 !!$ enddo 00197 !!$ enddo 00198 !!$ ! ERROR in type description!!! waves%global_param%f_assumption(nion+1) 00199 !!$ 00200 !!$ !print *, "itm_waves: global_param: nnphi", waves%global_param%nntor 00201 !!$ !print *, "itm_waves: global_param: nphi", waves%global_param%ntor 00202 !!$ !print *, "itm_waves: global_param: ptot", waves%global_param%power_tot 00203 !!$ !print *, "itm_waves: global_param: j", waves%global_param%cur_tor_ntor 00204 !!$ !print *, "itm_waves: global_param: freq", waves%global_param%frequency 00205 !!$ 00206 !!$ 00207 !!$ 00208 !!$ !---- waves%fullwave%local ----! 00209 !!$ allocate( waves%fullwave%local%e_plus(nfreq, max_nnphi, max_npsi,max_ntheta) ) 00210 !!$ allocate( waves%fullwave%local%e_minus(nfreq, max_nnphi, max_npsi,max_ntheta) ) 00211 !!$ allocate( waves%fullwave%local%e_para(nfreq, max_nnphi, max_npsi,max_ntheta) ) 00212 !!$ 00213 !!$ allocate( waves%fullwave%local%e_plus_ph(nfreq, max_nnphi, max_npsi,max_ntheta) ) 00214 !!$ allocate( waves%fullwave%local%e_minus_ph(nfreq, max_nnphi, max_npsi,max_ntheta) ) 00215 !!$ allocate( waves%fullwave%local%e_para_ph(nfreq, max_nnphi, max_npsi,max_ntheta) ) 00216 !!$ 00217 !!$ !---- waves%fullwave%local ----! 00218 !!$ do jfreq = 1, nfreq 00219 !!$ do jnphi = 1, waves%global_param%nntor(jfreq) 00220 !!$ do jpsi = 1, max_npsi 00221 !!$ do jtheta = 1, max_ntheta 00222 !!$ waves%fullwave%local%e_plus( jfreq, jnphi, jpsi,jtheta) = 0.5d0 00223 !!$ waves%fullwave%local%e_minus( jfreq, jnphi, jpsi,jtheta) = 0.5d0 00224 !!$ waves%fullwave%local%e_para( jfreq, jnphi, jpsi,jtheta) = 0d0 00225 !!$ waves%fullwave%local%e_plus_ph( jfreq, jnphi, jpsi,jtheta) = 0d0 00226 !!$ waves%fullwave%local%e_minus_ph(jfreq, jnphi, jpsi,jtheta) = 0d0 00227 !!$ waves%fullwave%local%e_para_ph( jfreq, jnphi, jpsi,jtheta) = 0d0 00228 !!$ enddo 00229 !!$ enddo 00230 !!$ enddo 00231 !!$ enddo 00232 !!$ 00233 !!$ !print *, "e+ = ",waves%fullwave%local%e_plus 00234 !!$ !print *, "e- = ",waves%fullwave%local%e_minus 00235 !!$ !print *, "ph e+ = ",waves%fullwave%local%e_plus_ph 00236 !!$ !print *, "ph e- = ",waves%fullwave%local%e_minus_ph 00237 !!$ 00238 !!$ 00239 !!$ end subroutine construct_dummy_euitm_waves 00240 00241 subroutine destruct_dummy_euitm_waves(waves) 00242 00243 type(type_waves), intent(inout) :: waves 00244 00245 ! Local 00246 integer :: jfreq 00247 integer :: nfreq 00248 00249 nfreq = size(waves%coherentwave) 00250 00251 !---- waves%global_param ----! 00252 do jfreq = 1, nfreq 00253 deallocate(waves%coherentwave(jfreq)%global_param%ntor) 00254 deallocate(waves%coherentwave(jfreq)%global_param%p_frac_ntor) 00255 deallocate(waves%coherentwave(jfreq)%global_param%pow_i) 00256 deallocate(waves%coherentwave(jfreq)%global_param%pow_ntor_i) 00257 deallocate(waves%coherentwave(jfreq)%global_param%pow_ntor_e) 00258 deallocate(waves%coherentwave(jfreq)%global_param%cur_tor_ntor) 00259 enddo 00260 00261 00262 !---- waves%fullwave%local ----! 00263 do jfreq = 1, nfreq 00264 deallocate( waves%coherentwave(jfreq)%fullwave%local%e_plus) 00265 deallocate( waves%coherentwave(jfreq)%fullwave%local%e_minus) 00266 deallocate( waves%coherentwave(jfreq)%fullwave%local%e_para) 00267 00268 deallocate( waves%coherentwave(jfreq)%fullwave%local%e_plus_ph) 00269 deallocate( waves%coherentwave(jfreq)%fullwave%local%e_minus_ph) 00270 deallocate( waves%coherentwave(jfreq)%fullwave%local%e_para_ph) 00271 enddo 00272 00273 deallocate(waves%coherentwave) 00274 00275 00276 end subroutine destruct_dummy_euitm_waves 00277 00278 end module euitm_waves_interface