aboutsummaryrefslogtreecommitdiff
path: root/Source/AcceleratorLattice/AcceleratorLattice.H
blob: b8151a1a0491e41a0e55dfe739682b9ce579eeb6 (plain) (blame)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
/* Copyright 2022 David Grote
 *
 * This file is part of WarpX.
 *
 * License: BSD-3-Clause-LBNL
 */
#ifndef WARPX_ACCELERATORLATTICE_ACCELERATORLATTICE_H_
#define WARPX_ACCELERATORLATTICE_ACCELERATORLATTICE_H_

#include "LatticeElementFinder.H"
#include "LatticeElements/Drift.H"
#include "LatticeElements/HardEdgedQuadrupole.H"
#include "LatticeElements/HardEdgedPlasmaLens.H"

#include <memory>
#include <string>

/* \brief AcceleratorLattice contains the lattice elements
 */
class AcceleratorLattice
{
public:

    AcceleratorLattice ();

    /**
     * \brief Read in the MAD like lattice description
     *
     * @param[in] root_name base name for the input variables
     * @param[in] z_location starting z location for the lattice
     */
    void
    ReadLattice (std::string const & root_name, amrex::ParticleReal & z_location);

    /* Refinement level at which this accelerator lattice is defined */
    int m_level;

    /* Flags where any lattice elements have been specified by the user */
    bool m_lattice_defined = false;

    /**
     * \brief Initialize the element finder instance at the given level of refinement
     *
     * @param[in] lev the level of refinement
     * @param[in] ba the box array at the level of refinement
     * @param[in] dm the distribution map at the level of refinement
     */
    void InitElementFinder (int const lev, amrex::BoxArray const & ba, amrex::DistributionMapping const & dm);

    /**
     * \brief Update the element finder, needed when the simulation frame has moved relative to the lab frame
     *
     * @param[in] lev the level of refinement
     */
    void UpdateElementFinder (int const lev);

    /* The lattice element finder handles the lookup that finds the elements at the particle locations.
     * It should follow the same grid layout as the main grids.
     */
    std::unique_ptr<amrex::LayoutData<LatticeElementFinder>> m_element_finder;

    /**
     * \brief Return an instance of a lattice finder associated with the grid and that can be used
     * on the device (such as a GPU)
     *
     * @param[in] a_pti the grid where the finder is needed
     * @param[in] a_offset the particle offset since the finded needs information about the particles as well
     */
    LatticeElementFinderDevice GetFinderDeviceInstance (WarpXParIter const& a_pti, int const a_offset) const;

    /* All of the available lattice element types */
    Drift h_drift;
    HardEdgedQuadrupole h_quad;
    HardEdgedPlasmaLens h_plasmalens;

};

#endif // WARPX_ACCELERATORLATTICE_ACCELERATORLATTICE_H_