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 : #pragma once 10 : 11 : #include "ComplexTypes.h" 12 : #include "ElementUserObject.h" 13 : #include "FFTData.h" 14 : 15 : template <typename T> 16 : class FFTBufferBase; 17 : 18 : #define usingFFTBufferBaseMembers \ 19 : using FFTBufferBase<T>::_dim; \ 20 : using FFTBufferBase<T>::_grid; \ 21 : using FFTBufferBase<T>::_real_space_data; \ 22 : using FFTBufferBase<T>::_reciprocal_space_data; \ 23 : using FFTBufferBase<T>::_real_space_data_start; \ 24 : using FFTBufferBase<T>::_reciprocal_space_data_start; \ 25 : using FFTBufferBase<T>::_real_space_data_stride; \ 26 : using FFTBufferBase<T>::_reciprocal_space_data_stride; \ 27 : using FFTBufferBase<T>::_how_many 28 : 29 : /** 30 : * Generic FFT interleaved data buffer base class 31 : */ 32 : template <typename T> 33 : class FFTBufferBase : public ElementUserObject 34 : { 35 : using ComplexT = typename ComplexType<T>::type; 36 : 37 : public: 38 : static InputParameters validParams(); 39 : 40 : FFTBufferBase(const InputParameters & parameters); 41 : 42 42 : virtual void initialize() {} 43 : virtual void execute(); 44 42 : virtual void finalize() {} 45 0 : virtual void threadJoin(const UserObject &) {} 46 : 47 : ///@{ transforms with proper scaling guaranteed 48 : virtual void forward(); 49 : virtual void backward(); 50 : ///@} 51 : 52 : ///@{ transforms without proper scaling guaranteed 53 : virtual void forwardRaw() = 0; 54 : virtual void backwardRaw() = 0; 55 : ///@} 56 : 57 : ///@{ scaling required after respective transform 58 0 : virtual Real forwardScale() { return 1.0; } 59 93 : virtual Real backwardScale() { return forwardScale(); } 60 : ///@} 61 : 62 : ///@{ buffer access 63 18 : FFTData<T> & realSpace() { return _real_space_data; } 64 105 : FFTData<ComplexT> & reciprocalSpace() { return _reciprocal_space_data; } 65 18 : const FFTData<T> & realSpace() const { return _real_space_data; } 66 0 : const FFTData<ComplexT> & reciprocalSpace() const { return _reciprocal_space_data; } 67 : 68 : ///@{ real space data access by location 69 : const T & operator()(const Point & p) const; 70 : T & operator()(const Point & p); 71 : ///@} 72 : 73 : /// return the number of grid cells along each dimension without padding 74 15 : const std::vector<int> & grid() const { return _grid; } 75 : 76 : /// return the size of the box 77 0 : const Point & getBoxSize() const { return _box_size; } 78 : 79 : /// return the minimum coordinate corner 80 0 : const Point & getMinCorner() const { return _min_corner; } 81 : 82 : /// return the maximum coordinate corner 83 0 : const Point & getMaxCorner() const { return _max_corner; } 84 : 85 : /// return the buffer dimension 86 0 : const unsigned int & dim() const { return _dim; } 87 : 88 : /// return the buffer dimension 89 0 : const std::vector<Real> & kTable(unsigned int d) const 90 : { 91 : mooseAssert(d < _dim, "invalid kTable component"); 92 0 : return _k_table[d]; 93 : } 94 : 95 : protected: 96 : ///@{ mesh data 97 : MooseMesh & _mesh; 98 : unsigned int _dim; 99 : ///@} 100 : 101 : /// grid size for FFT (needs to be signed for FFTW) 102 : std::vector<int> _grid; 103 : 104 : ///@{ simulation box extents 105 : Point _min_corner; 106 : Point _max_corner; 107 : Point _box_size; 108 : ///@} 109 : 110 : /// FFT grid cell volume 111 : Real _cell_volume; 112 : 113 : ///@{ FFT data buffer and unpadded number of grid cells 114 : FFTData<T> _real_space_data; 115 : FFTData<ComplexT> _reciprocal_space_data; 116 : ///@} 117 : 118 : /// pointer to the start of the data 119 : Real * _real_space_data_start; 120 : Complex * _reciprocal_space_data_start; 121 : 122 : /// stride in units of double size 123 : std::ptrdiff_t _real_space_data_stride; 124 : std::ptrdiff_t _reciprocal_space_data_stride; 125 : 126 : /// optional moose sister variabe (to obtain IC from) 127 : std::vector<const VariableValue *> _moose_variable; 128 : 129 : /// pretabulated k-vector components 130 : std::vector<std::vector<Real>> _k_table; 131 : 132 : /// cache the howMany value 133 : const std::size_t _how_many; 134 : };