aboutsummaryrefslogtreecommitdiff
path: root/Source/FieldSolver/FiniteDifferenceSolver/EvolveG.cpp
blob: 0c971d5776b703e538d38c414ea6c0b9714b7eb5 (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
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
/* Copyright 2020
 *
 * This file is part of WarpX.
 *
 * License: BSD-3-Clause-LBNL
 */


#include "FiniteDifferenceSolver.H"
#ifndef WARPX_DIM_RZ
#   include "FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianYeeAlgorithm.H"
#   include "FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianCKCAlgorithm.H"
#   include "FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianNodalAlgorithm.H"
#else
#   include "FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CylindricalYeeAlgorithm.H"
#endif
#include "Utils/TextMsg.H"
#include "Utils/WarpXAlgorithmSelection.H"

#include <AMReX.H>
#include <AMReX_Array4.H>
#include <AMReX_Config.H>
#include <AMReX_Extension.H>
#include <AMReX_GpuContainers.H>
#include <AMReX_GpuControl.H>
#include <AMReX_GpuLaunch.H>
#include <AMReX_GpuQualifiers.H>
#include <AMReX_IndexType.H>
#include <AMReX_MFIter.H>
#include <AMReX_MultiFab.H>
#include <AMReX_REAL.H>

#include <AMReX_BaseFwd.H>

#include <array>
#include <memory>

using namespace amrex;

void FiniteDifferenceSolver::EvolveG (
    std::unique_ptr<amrex::MultiFab>& Gfield,
    std::array<std::unique_ptr<amrex::MultiFab>,3> const& Bfield,
    amrex::Real const dt)
{
#ifdef WARPX_DIM_RZ
    // TODO Implement G update equation in RZ geometry
    amrex::ignore_unused(Gfield, Bfield, dt);
#else
    // Select algorithm
    if (m_do_nodal)
    {
        EvolveGCartesian<CartesianNodalAlgorithm>(Gfield, Bfield, dt);
    }
    else if (m_fdtd_algo == ElectromagneticSolverAlgo::Yee)
    {
        EvolveGCartesian<CartesianYeeAlgorithm>(Gfield, Bfield, dt);
    }
    else if (m_fdtd_algo == ElectromagneticSolverAlgo::CKC)
    {
        EvolveGCartesian<CartesianCKCAlgorithm>(Gfield, Bfield, dt);
    }
    else
    {
        amrex::Abort(Utils::TextMsg::Err("EvolveG: unknown FDTD algorithm"));
    }
#endif
}

#ifndef WARPX_DIM_RZ

template<typename T_Algo>
void FiniteDifferenceSolver::EvolveGCartesian (
    std::unique_ptr<amrex::MultiFab>& Gfield,
    std::array<std::unique_ptr<amrex::MultiFab>,3> const& Bfield,
    amrex::Real const dt)
{

    amrex::Real constexpr c2 = PhysConst::c * PhysConst::c;

#ifdef AMREX_USE_OMP
#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
#endif

    // Loop over grids and over tiles within each grid
    for (amrex::MFIter mfi(*Gfield, TilingIfNotGPU()); mfi.isValid(); ++mfi)
    {
        // Extract field data for this grid/tile
        amrex::Array4<amrex::Real> const& G = Gfield->array(mfi);
        amrex::Array4<amrex::Real> const& Bx = Bfield[0]->array(mfi);
        amrex::Array4<amrex::Real> const& By = Bfield[1]->array(mfi);
        amrex::Array4<amrex::Real> const& Bz = Bfield[2]->array(mfi);

        // Extract stencil coefficients
        amrex::Real const* const AMREX_RESTRICT coefs_x = m_stencil_coefs_x.dataPtr();
        amrex::Real const* const AMREX_RESTRICT coefs_y = m_stencil_coefs_y.dataPtr();
        amrex::Real const* const AMREX_RESTRICT coefs_z = m_stencil_coefs_z.dataPtr();

        const int n_coefs_x = m_stencil_coefs_x.size();
        const int n_coefs_y = m_stencil_coefs_y.size();
        const int n_coefs_z = m_stencil_coefs_z.size();

        // Extract tilebox to loop over
        amrex::Box const& tf = mfi.tilebox(Gfield->ixType().toIntVect());

        // Loop over cells and update G
        amrex::ParallelFor(tf, [=] AMREX_GPU_DEVICE (int i, int j, int k)
        {
            G(i,j,k) += c2 * dt * (T_Algo::UpwardDx(Bx, coefs_x, n_coefs_x, i, j, k)
                                 + T_Algo::UpwardDy(By, coefs_y, n_coefs_y, i, j, k)
                                 + T_Algo::UpwardDz(Bz, coefs_z, n_coefs_z, i, j, k));
        });
    }
}

#endif