cold_spec Subroutine

public subroutine cold_spec()

Calculates cold D-alpha emission

Arguments

None

Calls

proc~~cold_spec~~CallsGraph proc~cold_spec cold_spec proc~ind2sub ind2sub proc~cold_spec->proc~ind2sub proc~store_photons store_photons proc~cold_spec->proc~store_photons proc~get_plasma get_plasma proc~cold_spec->proc~get_plasma interface~parallel_sum parallel_sum proc~cold_spec->interface~parallel_sum interface~randn randn proc~cold_spec->interface~randn proc~get_fields get_fields proc~store_photons->proc~get_fields proc~get_passive_grid_indices get_passive_grid_indices proc~store_photons->proc~get_passive_grid_indices proc~get_indices get_indices proc~store_photons->proc~get_indices proc~uvw_to_xyz uvw_to_xyz proc~store_photons->proc~uvw_to_xyz proc~spectrum spectrum proc~store_photons->proc~spectrum proc~get_position get_position proc~get_plasma->proc~get_position proc~in_plasma in_plasma proc~get_plasma->proc~in_plasma proc~get_plasma->proc~uvw_to_xyz proc~xyz_to_uvw xyz_to_uvw proc~get_plasma->proc~xyz_to_uvw proc~parallel_sum_i1 parallel_sum_i1 interface~parallel_sum->proc~parallel_sum_i1 proc~parallel_sum_i0 parallel_sum_i0 interface~parallel_sum->proc~parallel_sum_i0 proc~parallel_sum_i2 parallel_sum_i2 interface~parallel_sum->proc~parallel_sum_i2 proc~parallel_sum_d4 parallel_sum_d4 interface~parallel_sum->proc~parallel_sum_d4 proc~parallel_sum_d5 parallel_sum_d5 interface~parallel_sum->proc~parallel_sum_d5 proc~parallel_sum_d2 parallel_sum_d2 interface~parallel_sum->proc~parallel_sum_d2 proc~parallel_sum_d3 parallel_sum_d3 interface~parallel_sum->proc~parallel_sum_d3 proc~parallel_sum_d0 parallel_sum_d0 interface~parallel_sum->proc~parallel_sum_d0 proc~parallel_sum_d1 parallel_sum_d1 interface~parallel_sum->proc~parallel_sum_d1 proc~get_fields->proc~in_plasma proc~get_fields->proc~uvw_to_xyz proc~get_fields->proc~xyz_to_uvw proc~calc_perp_vectors calc_perp_vectors proc~get_fields->proc~calc_perp_vectors mpi_allreduce mpi_allreduce proc~parallel_sum_i1->mpi_allreduce proc~parallel_sum_i0->mpi_allreduce proc~cyl_to_xyz cyl_to_xyz proc~get_position->proc~cyl_to_xyz proc~parallel_sum_i2->mpi_allreduce proc~parallel_sum_d4->mpi_allreduce proc~in_plasma->proc~xyz_to_uvw proc~cyl_to_uvw cyl_to_uvw proc~in_plasma->proc~cyl_to_uvw interface~interpol_coeff interpol_coeff proc~in_plasma->interface~interpol_coeff proc~parallel_sum_d5->mpi_allreduce proc~parallel_sum_d2->mpi_allreduce proc~parallel_sum_d3->mpi_allreduce proc~parallel_sum_d0->mpi_allreduce proc~parallel_sum_d1->mpi_allreduce proc~cyl_interpol3d_coeff cyl_interpol3D_coeff interface~interpol_coeff->proc~cyl_interpol3d_coeff proc~interpol1d_coeff interpol1D_coeff interface~interpol_coeff->proc~interpol1d_coeff proc~interpol2d_coeff_arr interpol2D_coeff_arr interface~interpol_coeff->proc~interpol2d_coeff_arr proc~interpol2d_coeff interpol2D_coeff interface~interpol_coeff->proc~interpol2d_coeff proc~cyl_interpol3d_coeff_arr cyl_interpol3D_coeff_arr interface~interpol_coeff->proc~cyl_interpol3d_coeff_arr proc~interpol1d_coeff_arr interpol1D_coeff_arr interface~interpol_coeff->proc~interpol1d_coeff_arr proc~cyl_to_xyz->proc~uvw_to_xyz proc~cyl_to_xyz->proc~cyl_to_uvw proc~cyl_interpol3d_coeff->proc~interpol2d_coeff proc~interpol2d_coeff_arr->proc~interpol2d_coeff proc~cyl_interpol3d_coeff_arr->proc~cyl_interpol3d_coeff proc~cyl_interpol3d_coeff_arr->proc~interpol2d_coeff proc~interpol1d_coeff_arr->proc~interpol1d_coeff

Called by

proc~~cold_spec~~CalledByGraph proc~cold_spec cold_spec program~fidasim fidasim program~fidasim->proc~cold_spec

Contents

Source Code


Source Code

subroutine cold_spec
    !+ Calculates cold D-alpha emission
    integer :: ic, i, j, k, it, ncell
    real(Float64), dimension(3) :: ri, vhalo, random3
    integer,dimension(3) :: ind
    !! Determination of the CX probability
    type(LocalProfiles) :: plasma
    real(Float64) :: cold_photons
    integer :: n = 10000

    !$OMP PARALLEL DO schedule(dynamic,1) private(i,j,k,ic,ind, &
    !$OMP& cold_photons, ri, vhalo, random3, plasma)
    loop_over_cells: do ic = istart, spec_chords%ncell, istep
        call ind2sub(beam_grid%dims,spec_chords%cell(ic),ind)
        i = ind(1) ; j = ind(2) ; k = ind(3)
        ri = [beam_grid%xc(i), beam_grid%yc(j), beam_grid%zc(k)]

        call get_plasma(plasma, pos=ri)
        cold_photons = plasma%denn(3)*tables%einstein(2,3)
        if(cold_photons.le.0.0) cycle loop_over_cells

        do it=1, n
            !! Cold Spectra
            call randn(random3)
            vhalo = plasma%vrot + sqrt(plasma%ti*0.5/(v2_to_E_per_amu*inputs%ai))*random3
            call store_photons(ri, vhalo, cold_photons/n, spec%cold)
        enddo
    enddo loop_over_cells
    !$OMP END PARALLEL DO

#ifdef _MPI
    !! Combine Spectra
    call parallel_sum(spec%cold)
#endif

end subroutine cold_spec