コード例 #1
0
    def mesh_intersect_1(self, rays, mesh):
        """
        Find the intersection of rays with the mesh using the Möller–Trumbore
        algorithm.
        """
        profiler.start('mesh_intersect_1')
        O = rays['origin']
        D = rays['direction']

        m = rays['mask'].copy()
        X = np.full(D.shape, np.nan, dtype=np.float64)

        # Copying these makes the code easier to read,
        # but may increase memory usage for dense meshes.
        p0 = mesh['points'][mesh['faces'][..., 0], :]
        p1 = mesh['points'][mesh['faces'][..., 1], :]
        p2 = mesh['points'][mesh['faces'][..., 2], :]

        epsilon = 1e-15

        num_rays = len(m)
        hits = np.empty(num_rays, dtype=np.int)
        m_temp = np.empty(num_rays, dtype=bool)
        m_temp_2 = np.zeros(num_rays, dtype=bool)

        for ii in range(mesh['faces'].shape[0]):
            m_temp[:] = m
            edge1 = p1[ii, :] - p0[ii, :]
            edge2 = p2[ii, :] - p0[ii, :]
            h = np.cross(D[m_temp], edge2)
            f = np.einsum('i,ji->j', edge1, h, optimize=True)
            m_temp &= ~((f > -epsilon) & (f < epsilon))
            if not np.any(m_temp):
                continue

            f = 1.0 / f
            s = O - p0[ii, :]
            u = f * np.einsum('ij,ij->i', s, h, optimize=True)
            m_temp &= ~((u < 0.0) | (u > 1.0))
            if not np.any(m_temp):
                continue

            q = np.cross(s, edge1)
            v = f * np.einsum('ij,ij->i', D, q, optimize=True)
            m_temp &= ~((v < 0.0) | (u + v > 1.0))
            if not np.any(m_temp):
                continue

            t = f * np.einsum('i,ji->j', edge2, q, optimize=True)

            # Update overall hit array and hit mask.
            m_temp_2[m_temp] = m_temp[m_temp]
            hits[m_temp] = ii
            X[m_temp] = O[m_temp] + t[m_temp, None] * D[m_temp, :]

        # Update the mask not to include any rays that didn't hit the mesh.
        m &= m_temp_2
        profiler.stop('mesh_intersect_1')

        return X, m, hits
コード例 #2
0
    def bundle_generate(self, bundle_input):

        profiler.start("Bundle Input Generation")
        m = bundle_input['mask']

        # Attempt to generate the specified number of bundles, but throw out
        # bundles that our outside of the last closed flux surface.
        #
        # This loop was setup for VMEC. Here we could do this as a single
        # vectorized operation instead.
        rho = np.zeros(len(m[m]))
        for ii in range(len(m[m])):
            # convert from cartesian coordinates to normalized radial coordinate.
            profiler.start("Fluxspace from Realspace")
            rho[ii] = self.rho_from_car(bundle_input['origin'][m][ii, :])
            profiler.stop("Fluxspace from Realspace")

        # evaluate emissivity, temperature and velocity at each bundle location.
        bundle_input['temperature'][m] = self.get_temperature(rho) * self.param['temperature_scale']
        bundle_input['emissivity'][m] = self.get_emissivity(rho) * self.param['emissivity_scale']
        bundle_input['velocity'][m] = self.get_velocity(rho) * self.param['velocity_scale']

        fintest = np.isfinite(bundle_input['temperature'])

        m &= fintest

        profiler.stop("Bundle Input Generation")

        return bundle_input
コード例 #3
0
def raytrace(config):
    """
    Perform a series of ray tracing runs.

    Each run will rebuild all objects, reset the random seed and then
    perform the requested number of iterations.

    If the option 'save_images' is set, then images will be saved
    at the completion of each run. The saving of these run images
    is one reason to use this routine rather than just increasing
    the number of iterations: periodic outputs during long computations.

    Also see :func:`~xicsrt.xicsrt_multiprocessing.raytrace` for a
    multiprocessing version of this routine.
    """
    profiler.start('raytrace')
    
    # Update the default config with the user config.
    config = xicsrt_config.get_config(config)
    check_config(config)

    # Make local copies of some options.
    num_runs = config['general']['number_of_runs']
    random_seed = config['general']['random_seed']
    
    output_list = []
    
    for ii in range(num_runs):
        m_log.info('Starting run: {} of {}'.format(ii + 1, num_runs))
        config_run = deepcopy(config)
        config_run['general']['output_run_suffix'] = '{:04d}'.format(ii)

        # Make sure each run uses a unique random seed.
        if random_seed is not None:
            random_seed += ii
        config_run['general']['random_seed'] = random_seed
        
        iteration = raytrace_single(config_run, _internal=True)
        output_list.append(iteration)
        
    output = combine_raytrace(output_list)

    # Reset the configuration options that were unique to the individual runs.
    output['config']['general']['output_run_suffix'] = config['general']['output_run_suffix']
    output['config']['general']['random_seed'] = config['general']['random_seed']

    if config['general']['save_config']:
        xicsrt_io.save_config(output['config'])
    if config['general']['save_images']:
        xicsrt_io.save_images(output)
    if config['general']['save_results']:
        xicsrt_io.save_results(output)
    if config['general']['print_results']:
        print_raytrace(output)

    profiler.stop('raytrace')
    return output
コード例 #4
0
 def mesh_get_index(self, hits, faces):
     """
     Match faces to face indexes, with a loop over faces.
     """
     profiler.start('mesh_get_index')
     idx_hits = np.empty(hits.shape[0], dtype=np.int)
     for ii, ff in enumerate(faces):
         m_temp = np.all(np.equal(ff, hits), axis=1)
         idx_hits[m_temp] = ii
     profiler.stop('mesh_get_index')
     return idx_hits
コード例 #5
0
    def find_near_faces(self, X, mesh, mask):
        m = mask
        profiler.start('find_near_faces')
        idx = mesh['points_tree'].query(X[m])[1]

        faces_idx = np.zeros((8, len(m)), dtype=np.int)
        faces_mask = np.zeros((8, len(m)), dtype=np.bool)

        faces_idx[:, m] = mesh['p_faces_idx'][:, idx]
        faces_mask[:, m] = mesh['p_faces_mask'][:, idx]
        profiler.stop('find_near_faces')
        return faces_idx, faces_mask
コード例 #6
0
    def shape_fd(self, a, b, delta=None):
        profiler.start('finite difference')
        if delta is None: delta = 1e-8

        xyz, _ = self.torus(a, b)
        xyz1, _ = self.torus(a + delta, b)
        xyz2, _ = self.torus(a, b + delta)

        vec1 = xyz1 - xyz
        vec2 = xyz2 - xyz
        norm_fd = np.cross(vec1, vec2)
        norm_fd /= np.linalg.norm(norm_fd)
        profiler.stop('finite difference')
        return xyz, norm_fd
コード例 #7
0
def _sort_raytrace(input, max_lost=None):
    """
    Sort the rays into 'lost' and 'found' rays, then truncate
    the number of lost rays.
    """
    if max_lost is None:
        max_lost = 1000

    profiler.start('_sort_raytrace')

    output = dict()
    output['config'] = input['config']
    output['total'] = dict()
    output['total']['meta'] = dict()
    output['total']['image'] = dict()
    output['found'] = dict()
    output['found']['meta'] = dict()
    output['found']['history'] = dict()
    output['lost'] = dict()
    output['lost']['meta'] = dict()
    output['lost']['history'] = dict()

    output['total']['meta'] = input['meta']
    output['total']['image'] = input['image']

    if len(input['history']) > 0:
        key_opt_list = list(input['history'].keys())
        key_opt_last = key_opt_list[-1]

        w_found = np.flatnonzero(input['history'][key_opt_last]['mask'])
        w_lost = np.flatnonzero(np.invert(input['history'][key_opt_last]['mask']))

        # Save only a portion of the lost rays so that our lost history does
        # not become too large.
        max_lost = min(max_lost, len(w_lost))
        index_lost = np.arange(len(w_lost))
        np.random.shuffle(index_lost)
        w_lost = w_lost[index_lost[:max_lost]]

        for key_opt in key_opt_list:
            output['found']['history'][key_opt] = dict()
            output['lost']['history'][key_opt] = dict()

            for key_ray in input['history'][key_opt]:
                output['found']['history'][key_opt][key_ray] = input['history'][key_opt][key_ray][w_found]
                output['lost']['history'][key_opt][key_ray] = input['history'][key_opt][key_ray][w_lost]

    profiler.stop('_sort_raytrace')

    return output
コード例 #8
0
def _raytrace_iter(config, sources, optics):
    """ 
    Perform a single iteration of raytracing with the given sources and optics.
    The returned rays are unsorted.
    """
    profiler.start('_raytrace_iter')

    # Setup local names for a few config entries.
    # This is only to make the code below more readable.
    keep_meta    = config['general']['keep_meta']
    keep_images  = config['general']['keep_images']
    keep_history = config['general']['keep_history']

    m_log.debug('Generating rays')
    rays = sources.generate_rays(keep_history=keep_history)
    m_log.debug('Raytracing optics')
    rays = optics.trace(rays, keep_history=keep_history, keep_images=keep_images)

    # Combine sources and optics outputs.
    meta    = dict()
    image   = dict()
    history = dict()
    
    if keep_meta:
        for key in sources.meta:
            meta[key] = sources.meta[key]
        for key in optics.meta:
            meta[key] = optics.meta[key]
    
    if keep_images:
        for key in sources.image:
            image[key] = sources.image[key]
        for key in optics.image:
            image[key] = optics.image[key]    
    
    if keep_history:
        for key in sources.history:
            history[key] = sources.history[key]
        for key in optics.history:
            history[key] = optics.history[key]

    output = dict()
    output['config'] = config
    output['meta'] = meta
    output['image'] = image
    output['history'] = history
    
    profiler.stop('_raytrace_iter')
    return output
コード例 #9
0
 def find_point_faces(self, p_idx, faces, mask=None):
     """
     Find all of the the faces that include a given mesh point.
     """
     profiler.start('find_point_faces')
     if mask is None:
         mask = np.ones(p_idx.shape, dtype=np.bool)
     m = mask
     p_faces_idx = np.zeros((8, len(m)), dtype=np.int)
     p_faces_mask = np.zeros((8, len(m)), dtype=np.bool)
     for ii_p in p_idx:
         ii_f = np.nonzero(np.equal(faces, p_idx[ii_p]))[0]
         faces_num = len(ii_f)
         p_faces_idx[:faces_num, ii_p] = ii_f
         p_faces_mask[:faces_num, ii_p] = True
     profiler.stop('find_point_faces')
     return p_faces_idx, p_faces_mask
コード例 #10
0
    def mesh_interpolate(self, X, mesh, mask):
        profiler.start('mesh_interpolate')
        X[:, 2] = mesh['interp']['z'](X[:, 0], X[:, 1])

        normals = np.empty(X.shape)
        normals[:, 0] = mesh['interp']['normal_x'](X[:, 0], X[:, 1])
        normals[:, 1] = mesh['interp']['normal_y'](X[:, 0], X[:, 1])
        normals[:, 2] = mesh['interp']['normal_z'](X[:, 0], X[:, 1])

        profiler.start('normalize')
        normals = np.einsum('i,ij->ij',
                            1.0 / np.linalg.norm(normals, axis=1),
                            normals,
                            optimize=True)
        profiler.stop('normalize')

        profiler.stop('mesh_interpolate')
        return X, normals
コード例 #11
0
    def calculate_mesh(self, a, b):
        profiler.start('calculate_mesh')

        num_a = len(a)
        num_b = len(b)

        aa, bb = np.meshgrid(a, b, indexing='ij')

        xx = np.empty((num_a, num_b))
        yy = np.empty((num_a, num_b))
        zz = np.empty((num_a, num_b))

        normal_xx = np.empty((num_a, num_b))
        normal_yy = np.empty((num_a, num_b))
        normal_zz = np.empty((num_a, num_b))

        # ------------------------------------------------
        # Now calculate the xyz values at each grid point.
        for ii_a in range(num_a):
            for ii_b in range(num_b):

                a = aa[ii_a, ii_b]
                b = bb[ii_a, ii_b]

                # Temporary for development.
                if self.param['normal_method'] == 'analytic':
                    xyz, norm = self.shape(a, b)
                elif self.param['normal_method'] == 'fd':
                    xyz, norm = self.shape_fd(a, b)
                elif self.param['normal_method'] == 'jax':
                    xyz, norm = self.shape_jax(a, b)
                else:
                    raise Exception(f"normal_method {self.param['normal_method']} unknown.")

                xx[ii_a, ii_b] = xyz[0]
                yy[ii_a, ii_b] = xyz[1]
                zz[ii_a, ii_b] = xyz[2]
                normal_xx[ii_a, ii_b] = norm[0]
                normal_yy[ii_a, ii_b] = norm[1]
                normal_zz[ii_a, ii_b] = norm[2]

        profiler.stop('calculate_mesh')

        return xx, yy, zz, normal_xx, normal_yy, normal_zz
コード例 #12
0
    def generate_mesh(self, mesh_size=None):
        """
        This method creates the meshgrid for the crystal
        """
        profiler.start('generate_mesh')

        # --------------------------------
        # Setup the basic grid parameters.

        a_range = self.param['angle_major']
        b_range = self.param['angle_minor']

        num_a = mesh_size[0]
        num_b = mesh_size[1]

        self.log.debug(f'num_a, num_b: {num_a}, {num_b}, total: {num_a*num_b}')

        a = np.linspace(a_range[0], a_range[1], num_a)
        b = np.linspace(b_range[0], b_range[1], num_b)

        xx, yy, zz, normal_xx, normal_yy, normal_zz = \
            self.calculate_mesh(a, b)

        aa, bb = np.meshgrid(a, b, indexing='ij')
        angles_2d = np.stack((aa.flatten(), bb.flatten()), axis=0).T
        tri = Delaunay(angles_2d)

        # It's also possible to triangulate using the x,y coordinates.
        # This is not recommended unless there is some specific need.
        #
        # points_2d = np.stack((xx.flatten(), yy.flatten()), axis=0).T
        # tri = Delaunay(points_2d)

        points = np.stack((xx.flatten(), yy.flatten(), zz.flatten())).T
        normals = np.stack((normal_xx.flatten(), normal_yy.flatten(), normal_zz.flatten())).T

        faces = tri.simplices

        profiler.stop('generate_mesh')
        return points, faces, normals
コード例 #13
0
    def trace(self, rays, keep_meta=None, keep_history=None, keep_images=None):
        """
        Perform raytracing for each object in sequence.
        """
        if keep_meta is None: keep_meta = True
        if keep_history is None: keep_history = False
        if keep_images is None: keep_images = False

        profiler.start('Dispatcher: raytrace')

        for key, obj in self.objects.items():

            profiler.start('Dispatcher: trace_global')
            rays = obj.trace_global(rays)
            profiler.stop('Dispatcher: trace_global')

            if keep_meta:
                self.meta[key] = dict()
                self.meta[key]['num_out'] = np.sum(rays['mask'])

            if keep_history:
                self.history[key] = deepcopy(rays)

            if keep_images:
                profiler.start('Dispatcher: collect')
                self.image[key] = obj.make_image(rays)
                profiler.stop('Dispatcher: collect')

        profiler.stop('Dispatcher: raytrace')

        return rays
コード例 #14
0
    def mesh_initialize(self):
        """
        Pre-calculate a number of mesh properties that are
        needed in the other mesh methods.
        """
        profiler.start('mesh_initialize')

        dummy = self._mesh_precalc(self.param['mesh_points'],
                                   self.param['mesh_normals'],
                                   self.param['mesh_faces'])
        self.param['mesh'] = {}
        for key in dummy:
            self.param['mesh'][key] = dummy[key]

        if self.param['mesh_coarse_points'] is not None:
            dummy = self._mesh_precalc(self.param['mesh_coarse_points'],
                                       self.param['mesh_coarse_normals'],
                                       self.param['mesh_coarse_faces'])
            self.param['mesh_coarse'] = {}
            for key in dummy:
                self.param['mesh_coarse'][key] = dummy[key]

        profiler.stop('mesh_initialize')
コード例 #15
0
    def generate_rays(self):
        rays = RayArray()
        profiler.start('generate_rays')
        
        profiler.start('generate_origin')
        rays['origin'] = self.generate_origin()
        profiler.stop('generate_origin')

        profiler.start('generate_direction')
        rays['direction'] = self.generate_direction(rays['origin'])
        profiler.stop('generate_direction')

        profiler.start('generate_wavelength')
        rays['wavelength'] = self.generate_wavelength(rays['direction'])
        profiler.stop('generate_wavelength')

        profiler.start('generate_weight')
        rays['weight'] = self.generate_weight()
        profiler.stop('generate_weight')
        
        profiler.start('generate_mask')
        rays['mask'] = self.generate_mask()
        profiler.stop('generate_mask')
        
        profiler.start('filter_rays')
        rays = self.ray_filter(rays)
        profiler.stop('filter_rays')        
        
        profiler.stop('generate_rays')
        return rays
コード例 #16
0
def raytrace_single(config, _internal=False):
    """
    Perform a single raytrace run consisting of multiple iterations.

    If history is enabled, sort the rays into those that are detected and
    those that are lost (found and lost). The found ray history will be
    returned in full. The lost ray history will be truncated to allow
    analysis of lost ray pattern while still limiting memory usage.

    private keywords
    ================
    _internal : bool (False)
      Used when calling this function from `raytrace` as part of the execution
      of multiple runs. Controls how `history_max_lost` is handled along with
      how `save_config` and `save_results` are interpreted.
    """
    profiler.start('raytrace_single')

    # Update the default config with the user config.
    config = xicsrt_config.config_to_numpy(config)
    config = xicsrt_config.get_config(config)
    check_config(config)

    m_log.info('Seeding np.random with {}'.format(config['general']['random_seed']))
    np.random.seed(config['general']['random_seed'])

    num_iter = config['general']['number_of_iter']
    max_lost_iter = int(config['general']['history_max_lost']/num_iter)

    if _internal:
        max_lost_iter = max_lost_iter//config['general']['number_of_runs']

    # Setup the dispatchers.
    if 'filters' in config:
        m_log.debug("Creating filters")
        filters = Dispatcher(config, 'filters')
        filters.instantiate()
        filters.setup()
        filters.initialize()
        config['filters'] = filters.get_config()
    else:
        filters = None

    m_log.debug("Creating sources")
    sources = Dispatcher(config, 'sources')
    sources.instantiate()
    sources.apply_filters(filters)
    sources.setup()
    sources.check_param()
    sources.initialize()
    config['sources'] = sources.get_config()

    m_log.debug("Creating optics")
    optics = Dispatcher(config, 'optics')
    optics.instantiate()
    optics.apply_filters(filters)
    optics.setup()
    optics.check_param()
    optics.initialize()
    config['optics'] = optics.get_config()

    # Do the actual raytracing
    output_list = []
    for ii in range(num_iter):
        m_log.info('Starting iteration: {} of {}'.format(ii + 1, num_iter))

        single = _raytrace_iter(config, sources, optics)
        sorted = _sort_raytrace(single, max_lost=max_lost_iter)
        output_list.append(sorted)

    output = combine_raytrace(output_list)

    if _internal is False:
        if config['general']['print_results']:
            print_raytrace(output)
        if config['general']['save_config']:
            xicsrt_io.save_config(output['config'])
        if config['general']['save_results']:
            xicsrt_io.save_results(output)

    if config['general']['save_images']:
        xicsrt_io.save_images(output)

    profiler.stop('raytrace_single')
    # profiler.report()
    return output
コード例 #17
0
def combine_raytrace(input_list):
    """
    Produce a combined output from a list of raytrace outputs.
    """
    profiler.start('combine_raytrace')

    output = dict()
    output['config'] = input_list[0]['config']
    output['total'] = dict()
    output['total']['meta'] = dict()
    output['total']['image'] = dict()
    output['found'] = dict()
    output['found']['meta'] = dict()
    output['found']['history'] = dict()
    output['lost'] = dict()
    output['lost']['meta'] = dict()
    output['lost']['history'] = dict()

    num_iter = len(input_list)
    key_opt_list = list(input_list[0]['total']['meta'].keys())
    key_opt_last = key_opt_list[-1]

    # Combine the meta data.
    for key_opt in key_opt_list:
        output['total']['meta'][key_opt] = dict()
        key_meta_list = list(input_list[0]['total']['meta'][key_opt].keys())
        for key_meta in key_meta_list:
            output['total']['meta'][key_opt][key_meta] = 0
            for ii_iter in range(num_iter):
                output['total']['meta'][key_opt][key_meta] += input_list[ii_iter]['total']['meta'][key_opt][key_meta]

    # Combine the images.
    for key_opt in key_opt_list:
        if key_opt in input_list[0]['total']['image']:
            if input_list[0]['total']['image'][key_opt] is not None:
                output['total']['image'][key_opt] = np.zeros(input_list[0]['total']['image'][key_opt].shape)
                for ii_iter in range(num_iter):
                    output['total']['image'][key_opt] += input_list[ii_iter]['total']['image'][key_opt]
            else:
                output['total']['image'][key_opt] = None

    # Combine all the histories.
    if len(input_list[0]['found']['history']) > 0:
        final_num_found = 0
        final_num_lost = 0
        for ii_run in range(num_iter):
            final_num_found += len(input_list[ii_run]['found']['history'][key_opt_last]['mask'])
            final_num_lost += len(input_list[ii_run]['lost']['history'][key_opt_last]['mask'])

        rays_found_temp = RayArray()
        rays_found_temp.zeros(final_num_found)

        rays_lost_temp = RayArray()
        rays_lost_temp.zeros(final_num_lost)

        for key_opt in key_opt_list:
            output['found']['history'][key_opt] = rays_found_temp.copy()
            output['lost']['history'][key_opt] = rays_lost_temp.copy()

        index_found = 0
        index_lost = 0
        for ii_run in range(num_iter):
            num_found = len(input_list[ii_run]['found']['history'][key_opt_last]['mask'])
            num_lost = len(input_list[ii_run]['lost']['history'][key_opt_last]['mask'])

            for key_opt in key_opt_list:
                for key_ray in output['found']['history'][key_opt]:
                    output['found']['history'][key_opt][key_ray][index_found:index_found + num_found] = (
                        input_list[ii_run]['found']['history'][key_opt][key_ray][:])
                    output['lost']['history'][key_opt][key_ray][index_lost:index_lost + num_lost] = (
                        input_list[ii_run]['lost']['history'][key_opt][key_ray][:])

            index_found += num_found
            index_lost += num_lost

    profiler.stop('combine_raytrace')
    return output
コード例 #18
0
    def create_sources(self, bundle_input):
        """
        Generate rays from a list of bundles.

        bundle_input
          a list containing dictionaries containing the locations, emissivities,
          temperatures and velocitities and of all ray bundles to be emitted.
        """

        rays_list = []
        count_rays_in_bundle = []

        m = bundle_input['mask']

        # Check if the number of rays generated will exceed max ray limits.
        # This is only approximate since poisson statistics may be in use.

        predicted_rays = int(
            np.sum(
                bundle_input['emissivity'][m] * self.param['time_resolution'] *
                self.param['bundle_volume'] * bundle_input['solid_angle'][m] /
                (4 * np.pi) * self.param['volume'] /
                (self.param['bundle_count'] * self.param['bundle_volume'])))

        self.log.debug(f'Predicted rays: {predicted_rays:0.2e}')
        if predicted_rays > self.param['max_rays']:
            raise ValueError(
                f"Current settings will produce too many rays ({predicted_rays:0.2e}). "
                f"Please reduce integration time or adjust other parameters.")

        # Bundle generation loop
        for ii in range(self.param['bundle_count']):
            if not bundle_input['mask'][ii]:
                continue
            profiler.start("Ray Bundle Generation")
            source_config = dict()

            # Specially dependent parameters
            source_config['origin'] = bundle_input['origin'][ii]
            source_config['temperature'] = bundle_input['temperature'][ii]
            source_config['velocity'] = bundle_input['velocity'][ii]
            source_config['spread'] = bundle_input['spread'][ii]

            # Calculate the total number of photons to launch from this bundle
            # volume. Since the source can use poisson statistics, this should
            # be of floating point type.
            intensity = (bundle_input['emissivity'][ii] *
                         self.param['time_resolution'] *
                         self.param['bundle_volume'] *
                         bundle_input['solid_angle'][ii] / (4 * np.pi))

            # Scale the number of photons based on the number of bundles.
            #
            # Ultimately we allow bundle_volume and bundle_count to be
            # independent, which means that a bundle representing a volume in
            # the plasma can be launched from virtual volume of a different
            # size.
            #
            # In order to allow this while maintaining overall photon statistics
            # from the plasma, we normalize the intensity so that each bundle
            # represents a volume of plasma_volume/bundle_count.
            #
            # In doing so bundle_volume cancels out, but I am leaving the
            # calculation separate for clarity.
            intensity *= self.param['volume'] / (self.param['bundle_count'] *
                                                 self.param['bundle_volume'])

            source_config['intensity'] = intensity

            # constants
            source_config['xsize'] = self.param['voxel_size']
            source_config['ysize'] = self.param['voxel_size']
            source_config['zsize'] = self.param['voxel_size']
            source_config['zaxis'] = self.param['zaxis']
            source_config['xaxis'] = self.param['xaxis']
            source_config['target'] = self.param['target']
            source_config['mass_number'] = self.param['mass_number']
            source_config['wavelength_dist'] = self.param['wavelength_dist']
            source_config['wavelength'] = self.param['wavelength']
            source_config['wavelength_range'] = self.param['wavelength_range']
            source_config['linewidth'] = self.param['linewidth']
            source_config['angular_dist'] = self.param['angular_dist']
            source_config['use_poisson'] = self.param['use_poisson']

            #create ray bundle sources and generate bundled rays
            source = XicsrtSourceFocused(source_config)
            bundled_rays = source.generate_rays()

            rays_list.append(bundled_rays)
            count_rays_in_bundle.append(len(bundled_rays['mask']))

            profiler.stop("Ray Bundle Generation")

        profiler.start('Ray Bundle Collection')
        # append bundled rays together to form a single ray dictionary.
        # create the final ray dictionary
        total_rays = np.int(np.sum(count_rays_in_bundle))
        rays = dict()
        rays['origin'] = np.zeros((total_rays, 3), dtype=np.float64)
        rays['direction'] = np.zeros((total_rays, 3), dtype=np.float64)
        rays['wavelength'] = np.zeros((total_rays), dtype=np.float64)
        rays['weight'] = np.zeros((total_rays), dtype=np.float64)
        rays['mask'] = np.ones((total_rays), dtype=np.bool)

        index = 0
        for ii, num_rays in enumerate(count_rays_in_bundle):
            rays['origin'][index:index + num_rays] = rays_list[ii]['origin']
            rays['direction'][index:index +
                              num_rays] = rays_list[ii]['direction']
            rays['wavelength'][index:index +
                               num_rays] = rays_list[ii]['wavelength']
            rays['weight'][index:index + num_rays] = rays_list[ii]['weight']
            rays['mask'][index:index + num_rays] = rays_list[ii]['mask']
            index += num_rays
        profiler.stop('Ray Bundle Collection')

        if len(rays['mask']) == 0:
            raise ValueError(
                'No rays generated. Check plasma input parameters')

        self.log.debug('Bundles Generated:       {:0.4e}'.format(len(m[m])))
        self.log.debug('Rays per bundle, mean:   {:0.0f}'.format(
            np.mean(count_rays_in_bundle)))
        self.log.debug('Rays per bundle, median: {:0.0f}'.format(
            np.median(count_rays_in_bundle)))
        self.log.debug('Rays per bundle, max:    {:0d}'.format(
            np.max(count_rays_in_bundle)))
        self.log.debug('Rays per bundle, min:    {:0d}'.format(
            np.min(count_rays_in_bundle)))

        return rays
コード例 #19
0
    def mesh_intersect_2(
        self,
        rays,
        mesh,
        mask,
        faces_idx,
        faces_mask,
    ):
        """
        Check for ray intersection with a list of mesh faces.

        Programming Notes
        -----------------
        Because of the mesh indexing, the arrays here have different
        dimensions than in mesh_intersect_1, and need a different
        vectorization.

        At the moment I am using an less efficient mesh intersection
        method. This should be updated to use the same method as
        mesh_intersect_1, but with the proper vectorization.
        """

        profiler.start('mesh_intersect_2')

        O = rays['origin']
        D = rays['direction']

        m = mask.copy()
        X = np.full(D.shape, np.nan, dtype=np.float64)

        num_rays = len(m)
        hits = np.empty(num_rays, dtype=np.int)
        epsilon = 1e-15

        # Copying these makes the code easier to read,
        # but may increase memory usage for dense meshes.
        faces = mesh['faces'][faces_idx]
        p0 = mesh['points'][faces[..., 0], :]
        p1 = mesh['points'][faces[..., 1], :]
        p2 = mesh['points'][faces[..., 2], :]
        n = mesh['faces_normal'][faces_idx]

        # distance = np.dot((p0 - O), n) / np.dot(D, n)
        t0 = p0[:, :, :] - O[None, :, :]
        t1 = np.einsum('ijk, ijk -> ij', t0, n, optimize=True)
        t2 = np.einsum('jk, ijk -> ij', D, n, optimize=True)
        dist = t1 / t2
        t3 = np.einsum('jk,ij -> ijk', D, dist, optimize=True)
        intersect = t3 + O

        # Pre-calculate some vectors and do the calculation in a
        # single step in an attempt to optimize this calculation.
        a = intersect - p0
        b = intersect - p1
        c = intersect - p2

        diff = (np.linalg.norm(np.cross(b, c), axis=2) +
                np.linalg.norm(np.cross(c, a), axis=2) +
                np.linalg.norm(np.cross(a, b), axis=2) -
                np.linalg.norm(np.cross((p0 - p1), (p0 - p2)), axis=2))

        # .. ToDo:
        #    For now hard code the floating point tolerance.
        #    A better way of handling floating point errors is needed.
        test = (diff < 1e-10) & (dist >= 0) & faces_mask
        m &= np.any(test, axis=0)

        # idx_hits tells us which of the 8 faces had a hit.
        idx_hits = np.argmax(test[:, m], axis=0)
        # Now index the faces_idx to git the actual face number.
        hits[m] = faces_idx[idx_hits, m]
        X[m] = intersect[idx_hits, m, :]

        profiler.stop('mesh_intersect_2')

        return X, m, hits
コード例 #20
0
    def _mesh_precalc(self, points, normals, faces):
        profiler.start('_mesh_precalc')

        output = {}
        output['faces'] = faces
        output['points'] = points
        output['normals'] = normals

        if faces is None:
            tri = Delaunay(points[:, 0:2])
            faces = tri.simplices
            output['faces'] = faces

        if self.param['mesh_interpolate']:
            # Create a set of interpolators.
            # For now create these in 2D using the x and y locations.
            # This will not make sense for all geometries. In principal
            # the interpolation could be done in 3D, but this needs more
            # investigation and will ultimately be less accurate.
            #
            # It is recommended that mesh optics be built using a local
            # coordinate system that makes the x and y coordinates sensible
            # for 2d interpolation.
            profiler.start('Create Interpolators')
            interp = {}
            output['interp'] = interp
            interp['z'] = Interpolator(points[:, 0:2], points[:, 2].flatten())
            interp['normal_x'] = Interpolator(points[:, 0:2],
                                              normals[:, 0].flatten())
            interp['normal_y'] = Interpolator(points[:, 0:2],
                                              normals[:, 1].flatten())
            interp['normal_z'] = Interpolator(points[:, 0:2],
                                              normals[:, 2].flatten())
            profiler.stop('Create Interpolators')

        # Copying these makes the code easier to read,
        # but may increase memory usage for dense meshes.
        p0 = points[faces[..., 0], :]
        p1 = points[faces[..., 1], :]
        p2 = points[faces[..., 2], :]

        # Calculate the normals at each face.
        faces_center = np.mean(np.array([p0, p1, p2]), 0)
        faces_normal = np.cross((p0 - p1), (p2 - p1))
        faces_normal /= np.linalg.norm(faces_normal, axis=1)[:, None]
        output['faces_center'] = faces_center
        output['faces_normal'] = faces_normal

        # Generate a tree for the points.
        points_tree = cKDTree(points)
        output['points_tree'] = points_tree

        # Build up a lookup table for the faces around each point.
        # This is currently slow for large arrays.
        points_idx = np.arange(len(points))
        p_faces_idx, p_faces_mask = \
            self.find_point_faces(points_idx, faces)
        output['p_faces_idx'] = p_faces_idx
        output['p_faces_mask'] = p_faces_mask

        # centers_tree = cKDTree(points)
        # output['centers_tree'] = centers_tree

        profiler.stop('_mesh_precalc')
        return output