/* Copyright 2019-2020 Andrew Myers, Maxence Thevenet, Remi Lehe * Revathi Jambunathan * * This file is part of WarpX. * * License: BSD-3-Clause-LBNL */ #include "Utils/WarpXConst.H" #include "SpectralKSpace.H" #include using namespace amrex; using namespace Gpu; /* \brief Initialize k space object. * * \param realspace_ba Box array that corresponds to the decomposition * of the fields in real space (cell-centered ; includes guard cells) * \param dm Indicates which MPI proc owns which box, in realspace_ba. * \param realspace_dx Cell size of the grid in real space */ SpectralKSpace::SpectralKSpace( const BoxArray& realspace_ba, const DistributionMapping& dm, const RealVect realspace_dx ) : dx(realspace_dx) // Store the cell size as member `dx` { AMREX_ALWAYS_ASSERT_WITH_MESSAGE( realspace_ba.ixType()==IndexType::TheCellType(), "SpectralKSpace expects a cell-centered box."); // Create the box array that corresponds to spectral space BoxList spectral_bl; // Create empty box list // Loop over boxes and fill the box list for (int i=0; i < realspace_ba.size(); i++ ) { // For local FFTs, boxes in spectral space start at 0 in // each direction and have the same number of points as the // (cell-centered) real space box // TODO: this will be different for the hybrid FFT scheme Box realspace_bx = realspace_ba[i]; IntVect fft_size = realspace_bx.length(); // Because the spectral solver uses real-to-complex FFTs, we only // need the positive k values along the fastest axis // (first axis for AMReX Fortran-order arrays) in spectral space. // This effectively reduces the size of the spectral space by half // see e.g. the FFTW documentation for real-to-complex FFTs IntVect spectral_bx_size = fft_size; spectral_bx_size[0] = fft_size[0]/2 + 1; // Define the corresponding box Box spectral_bx = Box( IntVect::TheZeroVector(), spectral_bx_size - IntVect::TheUnitVector() ); spectral_bl.push_back( spectral_bx ); } spectralspace_ba.define( spectral_bl ); // Allocate the components of the k vector: kx, ky (only in 3D), kz bool only_positive_k; for (int i_dim=0; i_dim& k = k_comp[mfi]; // Allocate k to the right size int N = bx.length( i_dim ); k.resize( N ); // Fill the k vector IntVect fft_size = realspace_ba[mfi].length(); const Real dk = 2*MathConst::pi/(fft_size[i_dim]*dx[i_dim]); AMREX_ALWAYS_ASSERT_WITH_MESSAGE( bx.smallEnd(i_dim) == 0, "Expected box to start at 0, in spectral space."); AMREX_ALWAYS_ASSERT_WITH_MESSAGE( bx.bigEnd(i_dim) == N-1, "Expected different box end index in spectral space."); if (only_positive_k){ // Fill the full axis with positive k values // (typically: first axis, in a real-to-complex FFT) for (int i=0; i& k = k_vec[i_dim][mfi]; ManagedVector& shift = shift_factor[mfi]; // Allocate shift coefficients shift.resize( k.size() ); // Fill the shift coefficients Real sign = 0; switch (shift_type){ case ShiftType::TransformFromCellCentered: sign = -1.; break; case ShiftType::TransformToCellCentered: sign = 1.; } const Complex I{0,1}; for (int i=0; i stencil_coef = getFonbergStencilCoefficients(n_order, nodal); // Loop over boxes and allocate the corresponding ManagedVector // for each box owned by the local MPI proc for ( MFIter mfi(spectralspace_ba, dm); mfi.isValid(); ++mfi ){ Real delta_x = dx[i_dim]; const ManagedVector& k = k_vec[i_dim][mfi]; ManagedVector& modified_k = modified_k_comp[mfi]; // Allocate modified_k to the same size as k modified_k.resize( k.size() ); // Fill the modified k vector for (int i=0; i getFonbergStencilCoefficients( const int n_order, const bool nodal ) { AMREX_ALWAYS_ASSERT_WITH_MESSAGE( n_order%2 == 0, "n_order should be even."); const int m = n_order/2; Vector coefs; coefs.resize( m+1 ); // Note: there are closed-form formula for these coefficients, // but they result in an overflow when evaluated numerically. // One way to avoid the overflow is to calculate the coefficients // by recurrence. // Coefficients for nodal (a.k.a. centered) finite-difference if (nodal == true) { coefs[0] = -2.; // First coefficient for (int n=1; n