From 960d1eae5714e785095e548fafe63ca94f687704 Mon Sep 17 00:00:00 2001 From: Daniel Schwen Date: Fri, 24 Jan 2020 14:51:36 -0700 Subject: [PATCH] Compute gradient (#401) --- examples/spectral.i | 54 ++++++++++++- .../executioners/SpectralExecutionerBase.h | 4 + include/userobjects/FFTBufferBase.h | 32 +++++++- src/executioners/SpectralExecutionerBase.C | 78 ++++++++++++++++++- src/userobjects/FFTBufferBase.C | 71 ++++++++++++++++- 5 files changed, 231 insertions(+), 8 deletions(-) diff --git a/examples/spectral.i b/examples/spectral.i index 437c565f..4ff3f251 100644 --- a/examples/spectral.i +++ b/examples/spectral.i @@ -45,6 +45,29 @@ function = 'cos(x/100*2*pi*4)*cos(y/100*2*pi*2)' [../] [../] + + + [./u_aux] + order = CONSTANT + family = MONOMIAL + [./InitialCondition] + type = SmoothCircleIC + x1 = 50 + y1 = 50 + radius = 30 + int_width = 20 + invalue = 1 + outvalue = 0 + [../] + [../] + [./grad_u0_aux] + order = CONSTANT + family = MONOMIAL + [../] + [./grad_u1_aux] + order = CONSTANT + family = MONOMIAL + [../] [] [Materials] @@ -66,8 +89,13 @@ moose_variable = 'R0_aux R1_aux R2_aux' [../] - # Solver - # ... + [./u] + type = RealFFTWBuffer + moose_variable = u_aux + [../] + [./grad_u] + type = RealVectorValueFFTWBuffer + [../] [] [AuxKernels] @@ -99,6 +127,28 @@ component = 2 execute_on = FINAL [../] + + [./u_aux] + type = FFTBufferAux + variable = u_aux + fft_buffer = u + execute_on = FINAL + [../] + + [./grad_u0_aux] + type = FFTBufferAux + variable = grad_u0_aux + fft_buffer = grad_u + component = 0 + execute_on = FINAL + [../] + [./grad_u1_aux] + type = FFTBufferAux + variable = grad_u1_aux + fft_buffer = grad_u + component = 1 + execute_on = FINAL + [../] [] [Executioner] diff --git a/include/executioners/SpectralExecutionerBase.h b/include/executioners/SpectralExecutionerBase.h index c690ad5d..3488383a 100644 --- a/include/executioners/SpectralExecutionerBase.h +++ b/include/executioners/SpectralExecutionerBase.h @@ -33,9 +33,13 @@ class SpectralExecutionerBase : public Executioner virtual bool lastSolveConverged() const override { return true; } protected: + /// obtain a non-const reference to an FFT buffer template FFTBufferBase & getFFTBuffer(const std::string & name); + /// multiply a scalar buffer by its corresponding k-vector field int a vector buffer + void kVectorMultiply(const FFTBufferBase & in, FFTBufferBase & out) const; + Real _system_time; int & _time_step; Real & _time; diff --git a/include/userobjects/FFTBufferBase.h b/include/userobjects/FFTBufferBase.h index 80cbb0fc..c7e33aab 100644 --- a/include/userobjects/FFTBufferBase.h +++ b/include/userobjects/FFTBufferBase.h @@ -53,6 +53,31 @@ class FFTBufferBase : public ElementUserObject T & operator()(const Point & p); ///@} + ///@{ convenience math operators + FFTBufferBase & operator+=(FFTBufferBase const & rhs); + FFTBufferBase & operator-=(FFTBufferBase const & rhs); + FFTBufferBase & operator*=(FFTBufferBase const & rhs); + FFTBufferBase & operator/=(FFTBufferBase const & rhs); + FFTBufferBase & operator*=(Real rhs); + FFTBufferBase & operator/=(Real rhs); + ///@} + + /// return the number of grid cells along each dimension without padding + const std::vector & grid() const { return _grid; } + + /// return the number of proper grid cells without padding + const std::size_t & size() const { return _grid_size; } + + /// return the buffer dimension + const unsigned int & dim() const { return _dim; } + + /// return the buffer dimension + const std::vector & kTable(unsigned int d) const + { + mooseAssert(d < _dim, "invalid kTable component"); + return _k_table[d]; + } + protected: /// get the addres of the first data element of the ith object in the bufefr Real * start(std::size_t i); @@ -77,9 +102,9 @@ class FFTBufferBase : public ElementUserObject /// FFT grid cell volume Real _cell_volume; - ///@{ FFT data buffer + ///@{ FFT data buffer and unpadded number of grid cells std::vector _buffer; - std::size_t _buffer_size; + std::size_t _grid_size; ///@} /// pointer to the start of the data @@ -91,6 +116,9 @@ class FFTBufferBase : public ElementUserObject /// optional moose sister variabe (to obtain IC from) std::vector _moose_variable; + /// pretabulated k-vector components + std::vector> _k_table; + /// cache the howMany value const std::size_t _how_many; }; diff --git a/src/executioners/SpectralExecutionerBase.C b/src/executioners/SpectralExecutionerBase.C index 13ad4c78..e7b0bdbb 100644 --- a/src/executioners/SpectralExecutionerBase.C +++ b/src/executioners/SpectralExecutionerBase.C @@ -56,23 +56,99 @@ SpectralExecutionerBase::execute() _fe_problem.outputStep(EXEC_INITIAL); _fe_problem.advanceState(); + // back and forth test auto & c = getFFTBuffer("c"); c.forward(); - auto & R = getFFTBuffer("R"); R.forward(); + // gradient test + auto & u = getFFTBuffer("u"); + u.forward(); + auto & grad_u = getFFTBuffer("grad_u"); + kVectorMultiply(u, grad_u); + _time_step = 1; _fe_problem.execute(EXEC_FINAL); _time = _time_step; _fe_problem.outputStep(EXEC_FINAL); _fe_problem.advanceState(); + // back and forth test c.backward(); R.backward(); + R /= 10000.0; + c /= 10000.0; + + // gradient test + u.backward(); + grad_u.backward(); + u /= 10000.0; + grad_u /= 100.0; _time_step = 2; _fe_problem.execute(EXEC_FINAL); _time = _time_step; _fe_problem.outputStep(EXEC_FINAL); } + +void +SpectralExecutionerBase::kVectorMultiply(const FFTBufferBase & in, + FFTBufferBase & out) const +{ + mooseAssert(in.size() == out.size(), "Buffer sizes must be equal"); + mooseAssert(in.dim() == out.dim(), "Buffer dimensions must be equal"); + + const auto & grid = in.grid(); + switch (in.dim()) + { + case 1: + { + const int ni = grid[0]; + for (int i = 0; i < ni; ++i) + { + out[i](0) = in[i] * i; + out[i](1) = 0.0; + out[i](2) = 0.0; + } + return; + } + + case 2: + { + std::size_t index = 0; + const int ni = grid[0]; + const int nj = grid[1]; + for (int i = 0; i < ni; ++i) + for (int j = 0; j < nj; ++j) + { + out[index](0) = in[index] * i; + out[index](1) = in[index] * j; + out[index](2) = 0.0; + index++; + } + return; + } + + case 3: + { + std::size_t index = 0; + const int ni = grid[0]; + const int nj = grid[1]; + const int nk = grid[2]; + for (int i = 0; i < ni; ++i) + for (int j = 0; j < nj; ++j) + for (int k = 0; k < nk; ++k) + { + out[index](0) = in[index] * i; + out[index](1) = in[index] * j; + out[index](2) = in[index] * k; + index++; + } + return; + } + + default: + mooseError("Invalid buffer dimension"); + } +} diff --git a/src/userobjects/FFTBufferBase.C b/src/userobjects/FFTBufferBase.C index e3b1f47b..af984c19 100644 --- a/src/userobjects/FFTBufferBase.C +++ b/src/userobjects/FFTBufferBase.C @@ -44,8 +44,9 @@ FFTBufferBase::FFTBufferBase(const InputParameters & parameters) _mesh(_subproblem.mesh()), _dim(_mesh.dimension()), _cell_volume(1.0), - _buffer_size(1), + _grid_size(1), _moose_variable(coupledComponents("moose_variable")), + _k_table(_dim), _how_many(howMany()) { // make sure Real is double @@ -98,6 +99,7 @@ FFTBufferBase::FFTBufferBase(const InputParameters & parameters) } // get mesh extents and calculate space required and estimate spectrum bins + std::size_t buffer_size = 1; for (unsigned int i = 0; i < _dim; ++i) { _min_corner(i) = _mesh.getMinInDimension(i); @@ -105,10 +107,18 @@ FFTBufferBase::FFTBufferBase(const InputParameters & parameters) _box_size(i) = _max_corner(i) - _min_corner(i); _cell_volume *= _box_size(i) / _grid[i]; + // unpadded size (number of physical grid cells) + _grid_size *= _grid[i]; + // last direction needs to be padded for in-place transforms - _buffer_size *= (i == _dim - 1) ? ((_grid[i] >> 1) + 1) << 1 : _grid[i]; + buffer_size *= (i == _dim - 1) ? ((_grid[i] >> 1) + 1) << 1 : _grid[i]; + + // precompute kvector components for current direction + _k_table[i].resize(_grid[i]); + for (int j = 0; j < _grid[i]; ++j) + _k_table[i][j] = (j * 2 > _grid[i] ? Real(_grid[i] - j) : Real(j)) / _box_size(i); } - _buffer.resize(_buffer_size); + _buffer.resize(buffer_size); // compute stride and start pointer _start = reinterpret_cast(start(0)); @@ -153,6 +163,61 @@ FFTBufferBase::operator()(const Point & p) return _buffer[a]; } +template +FFTBufferBase & +FFTBufferBase::operator+=(FFTBufferBase const & rhs) +{ + for (std::size_t i = 0; i < _grid_size; ++i) + _buffer[i] += rhs[i]; + return *this; +} + +template +FFTBufferBase & +FFTBufferBase::operator-=(FFTBufferBase const & rhs) +{ + for (std::size_t i = 0; i < _grid_size; ++i) + _buffer[i] -= rhs[i]; + return *this; +} + +template +FFTBufferBase & +FFTBufferBase::operator*=(FFTBufferBase const & rhs) +{ + for (std::size_t i = 0; i < _grid_size; ++i) + _buffer[i] *= rhs[i]; + return *this; +} + +template +FFTBufferBase & +FFTBufferBase::operator/=(FFTBufferBase const & rhs) +{ + for (std::size_t i = 0; i < _grid_size; ++i) + _buffer[i] /= rhs[i]; + return *this; +} + +template +FFTBufferBase & +FFTBufferBase::operator*=(Real rhs) +{ + for (std::size_t i = 0; i < _grid_size; ++i) + _buffer[i] *= rhs; + return *this; +} + +template +FFTBufferBase & +FFTBufferBase::operator/=(Real rhs) +{ + const Real reciprocal = 1 / rhs; + for (std::size_t i = 0; i < _grid_size; ++i) + _buffer[i] *= reciprocal; + return *this; +} + template <> Real * FFTBufferBase::start(std::size_t i)