halo_spec Subroutine

public subroutine halo_spec()

Calculates halo emission from user supplied neutrals file

Arguments

None

Calls

proc~~halo_spec~~CallsGraph proc~halo_spec halo_spec proc~mc_halo mc_halo proc~halo_spec->proc~mc_halo proc~ind2sub ind2sub proc~halo_spec->proc~ind2sub proc~get_plasma get_plasma proc~halo_spec->proc~get_plasma proc~store_photons store_photons proc~halo_spec->proc~store_photons interface~randu randu proc~halo_spec->interface~randu proc~in_plasma in_plasma proc~halo_spec->proc~in_plasma interface~parallel_sum parallel_sum proc~halo_spec->interface~parallel_sum interface~randn randn proc~mc_halo->interface~randn proc~get_plasma->proc~in_plasma proc~get_position get_position proc~get_plasma->proc~get_position proc~uvw_to_xyz uvw_to_xyz proc~get_plasma->proc~uvw_to_xyz proc~xyz_to_uvw xyz_to_uvw proc~get_plasma->proc~xyz_to_uvw proc~get_passive_grid_indices get_passive_grid_indices proc~store_photons->proc~get_passive_grid_indices proc~spectrum spectrum proc~store_photons->proc~spectrum proc~get_indices get_indices proc~store_photons->proc~get_indices proc~get_fields get_fields proc~store_photons->proc~get_fields proc~store_photons->proc~uvw_to_xyz proc~cyl_to_uvw cyl_to_uvw proc~in_plasma->proc~cyl_to_uvw proc~in_plasma->proc~xyz_to_uvw interface~interpol_coeff interpol_coeff proc~in_plasma->interface~interpol_coeff proc~parallel_sum_d1 parallel_sum_d1 interface~parallel_sum->proc~parallel_sum_d1 proc~parallel_sum_i2 parallel_sum_i2 interface~parallel_sum->proc~parallel_sum_i2 proc~parallel_sum_i3 parallel_sum_i3 interface~parallel_sum->proc~parallel_sum_i3 proc~parallel_sum_d3 parallel_sum_d3 interface~parallel_sum->proc~parallel_sum_d3 proc~parallel_sum_d2 parallel_sum_d2 interface~parallel_sum->proc~parallel_sum_d2 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_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_d0 parallel_sum_d0 interface~parallel_sum->proc~parallel_sum_d0 proc~xyz_to_cyl xyz_to_cyl proc~get_passive_grid_indices->proc~xyz_to_cyl proc~cyl_to_xyz cyl_to_xyz proc~get_position->proc~cyl_to_xyz proc~parallel_sum_d1->proc~parallel_sum_d1 mpi_allreduce mpi_allreduce proc~parallel_sum_d1->mpi_allreduce proc~parallel_sum_i2->proc~parallel_sum_i1 proc~parallel_sum_i2->mpi_allreduce proc~parallel_sum_i3->proc~parallel_sum_i1 proc~parallel_sum_i3->mpi_allreduce proc~parallel_sum_d3->proc~parallel_sum_d1 proc~parallel_sum_d3->mpi_allreduce proc~get_fields->proc~in_plasma proc~get_fields->proc~get_position 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 proc~parallel_sum_d2->proc~parallel_sum_d1 proc~parallel_sum_d2->mpi_allreduce proc~parallel_sum_d4->proc~parallel_sum_d1 proc~parallel_sum_d4->mpi_allreduce proc~parallel_sum_d5->proc~parallel_sum_d1 proc~parallel_sum_d5->mpi_allreduce proc~parallel_sum_i1->proc~parallel_sum_i1 proc~parallel_sum_i1->mpi_allreduce proc~parallel_sum_i0->mpi_allreduce proc~parallel_sum_d0->mpi_allreduce proc~xyz_to_cyl->proc~xyz_to_uvw proc~uvw_to_cyl uvw_to_cyl proc~xyz_to_cyl->proc~uvw_to_cyl proc~cyl_to_xyz->proc~cyl_to_uvw proc~cyl_to_xyz->proc~uvw_to_xyz

Called by

proc~~halo_spec~~CalledByGraph proc~halo_spec halo_spec program~fidasim fidasim program~fidasim->proc~halo_spec

Contents

Source Code


Source Code

subroutine halo_spec
    !+ Calculates halo emission from user supplied neutrals file
    integer :: ic, i, j, k, it, is
    real(Float64), dimension(3) :: ri, vhalo, random3, rc
    integer,dimension(3) :: ind
    !! Determination of the CX probability
    type(LocalProfiles) :: plasma
    real(Float64) :: halo_photons, wght
    logical :: inp
    integer :: n = 10000

    !$OMP PARALLEL DO schedule(dynamic,1) private(i,j,k,ic,ind,is, &
    !$OMP& halo_photons, wght, rc, ri, inp, 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)

        halo_photons = neut%halo%dens(3,i,j,k)*tables%einstein(2,3)
        if(halo_photons.le.0.0) cycle loop_over_cells

        rc = [beam_grid%xc(i), beam_grid%yc(j), beam_grid%zc(k)]

        !Find a point in cell that is also in the plasma
        ri = rc
        call in_plasma(ri, inp)
        do while (.not.inp)
            call randu(random3)
            ri = rc + beam_grid%dr*(random3 - 0.5)
            call in_plasma(ri,inp)
        enddo

        call get_plasma(plasma, pos=ri)
        do is=1, n_thermal
            wght = plasma%deni(is)/sum(plasma%deni)
            do it=1, n
                !! Halo Spectra
                call mc_halo(plasma, thermal_mass(is), vhalo)
                call store_photons(ri, vhalo, thermal_lambda0(is), wght*halo_photons/n, &
                    &spec%halo(:,:,:,is),spec%halostokes(:,:,:,:,is))
            enddo
        enddo
    enddo loop_over_cells
    !$OMP END PARALLEL DO

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

end subroutine halo_spec