Line data Source code
1 : /**********************************************************************/ 2 : /* DO NOT MODIFY THIS HEADER */ 3 : /* MAGPIE - Mesoscale Atomistic Glue Program for Integrated Execution */ 4 : /* */ 5 : /* Copyright 2017 Battelle Energy Alliance, LLC */ 6 : /* ALL RIGHTS RESERVED */ 7 : /**********************************************************************/ 8 : 9 : #include "SpectralExecutionerBase.h" 10 : #include "InputParameters.h" 11 : 12 : // testing 13 : registerMooseObject("MagpieApp", SpectralExecutionerBase); 14 : 15 : InputParameters 16 18 : SpectralExecutionerBase::validParams() 17 : { 18 18 : InputParameters params = Executioner::validParams(); 19 18 : params.addClassDescription("Executioner for FFT simulations."); 20 36 : params.addParam<Real>("time", 0.0, "System time"); 21 18 : return params; 22 0 : } 23 : 24 9 : SpectralExecutionerBase::SpectralExecutionerBase(const InputParameters & parameters) 25 : : Executioner(parameters), 26 9 : _system_time(getParam<Real>("time")), 27 9 : _time_step(_fe_problem.timeStep()), 28 9 : _time(_fe_problem.time()), 29 9 : _final_timer(registerTimedSection("final", 1)), 30 18 : _fft_problem(dynamic_cast<FFTProblem *>(&_fe_problem)) 31 : { 32 : 33 9 : if (!_fft_problem) 34 0 : mooseError("Use Problem/type=FFTProblem with a spectral executioner"); 35 9 : } 36 : 37 : void 38 9 : SpectralExecutionerBase::init() 39 : { 40 9 : if (_app.isRecovering()) 41 : { 42 0 : _console << "\nCannot recover FFT solves!\nExiting...\n" << std::endl; 43 0 : return; 44 : } 45 : 46 : // checkIntegrity(); 47 9 : _fe_problem.execute(EXEC_PRE_MULTIAPP_SETUP); 48 9 : _fe_problem.initialSetup(); 49 : } 50 : 51 : void 52 0 : SpectralExecutionerBase::execute() 53 : { 54 0 : mooseError("SpectralExecutionerBase is an abstract base class"); 55 : } 56 : 57 : void 58 6 : SpectralExecutionerBase::kVectorMultiply(const FFTBufferBase<Real> & in_buffer, 59 : FFTBufferBase<RealVectorValue> & out_buffer) const 60 : { 61 : mooseAssert(in_buffer.dim() == out_buffer.dim(), "Buffer dimensions must be equal"); 62 : 63 : const FFTData<Complex> & in = in_buffer.reciprocalSpace(); 64 : FFTData<ComplexVectorValue> & out = out_buffer.reciprocalSpace(); 65 : mooseAssert(in.size() == out.size(), "Buffer sizes must be equal"); 66 : 67 : const auto & grid = in_buffer.grid(); 68 : const Complex I(0.0, 1.0); 69 6 : switch (in_buffer.dim()) 70 : { 71 : case 1: 72 : { 73 : const auto & ivec = in_buffer.kTable(0); 74 0 : const int ni = grid[0]; 75 0 : for (int i = 0; i * 2 <= ni; ++i) 76 0 : out[i](0) = in[i] * ivec[i] * I; 77 : return; 78 : } 79 : 80 6 : case 2: 81 : { 82 : std::size_t index = 0; 83 : const auto & ivec = in_buffer.kTable(0); 84 : const auto & jvec = in_buffer.kTable(1); 85 6 : const int ni = grid[0]; 86 6 : const int nj = (grid[1] >> 1) + 1; 87 606 : for (int i = 0; i < ni; ++i) 88 16200 : for (int j = 0; j < nj; ++j) 89 : { 90 15600 : out[index](0) = in[index] * ivec[i] * I; 91 15600 : out[index](1) = in[index] * jvec[j] * I; 92 15600 : index++; 93 : } 94 : return; 95 : } 96 : 97 0 : case 3: 98 : { 99 : std::size_t index = 0; 100 : const auto & ivec = in_buffer.kTable(0); 101 : const auto & jvec = in_buffer.kTable(1); 102 : const auto & kvec = in_buffer.kTable(2); 103 0 : const int ni = grid[0]; 104 0 : const int nj = grid[1]; 105 0 : const int nk = grid[2]; 106 0 : for (int i = 0; i < ni; ++i) 107 0 : for (int j = 0; j < nj; ++j) 108 0 : for (int k = 0; k * 2 <= nk; ++k) 109 : { 110 0 : out[index](0) = in[index] * ivec[i] * I; 111 0 : out[index](1) = in[index] * jvec[j] * I; 112 0 : out[index](2) = in[index] * kvec[k] * I; 113 0 : index++; 114 : } 115 : return; 116 : } 117 : 118 0 : default: 119 0 : mooseError("Invalid buffer dimension"); 120 : } 121 : }