1#ifndef RPG_FOURTH_ORDER_PARAMETER_TPP
2#define RPG_FOURTH_ORDER_PARAMETER_TPP
4#include "FourthOrderParameter.h"
6#include <rpg/fts/simulator/Simulator.h>
9#include <prdc/crystal/shiftToMinimum.h>
10#include <prdc/cuda/resources.h>
12#include <pscf/mesh/MeshIterator.h>
13#include <pscf/math/IntVec.h>
15#include <util/param/ParamComposite.h>
16#include <util/misc/FileMaster.h>
17#include <util/misc/ioUtil.h>
18#include <util/format/Int.h>
19#include <util/format/Dbl.h>
62 IntVec<D> const & dimensions = system().mesh().dimensions();
65 for (
int i = 0; i < D; ++i) {
67 kMeshDimensions_[i] = dimensions[i];
68 kSize_ *= dimensions[i];
70 kMeshDimensions_[i] = dimensions[i]/2 + 1;
71 kSize_ *= (dimensions[i]/2 + 1);
77 wc0_.allocate(dimensions);
78 wK_.allocate(dimensions);
79 prefactor_.allocate(kMeshDimensions_);
83 isInitialized_ =
true;
85 if (!isInitialized_) {
86 UTIL_THROW(
"Error: object is not initialized");
98 const int nMonomer = system().mixture().nMonomer();
101 if (!simulator().hasWc()){
102 simulator().computeWc();
105 const int meshSize = system().domain().mesh().size();
111 system().fft().forwardTransform(wc0_, wK_);
121 FourthOrderParameter_ = std::pow(FourthOrderParameter_, 0.25);
123 return FourthOrderParameter_;
129 if (simulator().hasRamp() && nSamplePerOutput() == 1) {
130 double chi= system().interaction().chi(0,1);
133 outputFile_ <<
Int(step);
134 outputFile_ <<
Dbl(chi);
135 outputFile_ <<
Dbl(value);
145 IntVec<D> meshDimensions = system().domain().mesh().dimensions();
146 UnitCell<D> const & unitCell = system().domain().unitCell();
148 for (
int i = 0; i < kSize_; ++i){
162 Gmin = shiftToMinimum(G, meshDimensions, unitCell);
163 GminList[itr.
rank()] = Gmin;
168 bool inverseFound =
false;
171 if (prefactor_h[itr.
rank()] == 0){
172 Gmin = GminList[itr.
rank()];
179 for (; !searchItr.
atEnd(); ++searchItr){
180 if (nGmin == GminList[searchItr.
rank()]){
181 prefactor_h[itr.
rank()] = 1.0/2.0;
182 prefactor_h[searchItr.
rank()] = 1.0/2.0;
187 if (inverseFound ==
false){
188 prefactor_h[itr.
rank()] = 1.0;
196 prefactor_ = prefactor_h;
Template for dynamic array stored in host CPU memory.
An IntVec<D, T> is a D-component vector of elements of integer type T.
Iterator over points in a Mesh<D>.
int rank() const
Get the rank of current element.
void begin()
Set iterator to the first point in the mesh.
bool atEnd() const
Is this the end (i.e., one past the last point)?
IntVec< D > position() const
Get current position in the grid, as integer vector.
Field of real double precision values on an FFT mesh.
void allocate(IntVec< D > const &meshDimensions)
Allocate the underlying C array for an FFT grid.
Base template for UnitCell<D> classes, D=1, 2 or 3.
Analyze averages and block averages of several real variables.
virtual void setup()
Setup before loop.
virtual void outputValue(int step, double value)
Output a sampled or block average value.
virtual ~FourthOrderParameter()
Destructor.
virtual double compute()
Compute and return the fourth order parameter.
void setClassName(const char *className)
Set class name string.
FourthOrderParameter(Simulator< D > &simulator, System< D > &system)
Constructor.
virtual void outputValue(int step, double value)
Output a sampled or block average value.
virtual void setup()
Setup before simulation loop.
void computePrefactor()
Compute prefactor for each Fourier wavevector.
Field theoretic simulator (base class).
Main class for calculations that represent one system.
Vec< D, T > & negate(const Vec< D, T > &v)
Return negative of vector v.
Dynamically allocatable contiguous array template.
void allocate(int capacity)
Allocate the underlying C array.
Wrapper for a double precision number, for formatted ostream output.
Wrapper for an int, for formatted ostream output.
File containing preprocessor macros for error handling.
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
#define UTIL_THROW(msg)
Macro for throwing an Exception, reporting function, file and line number.
cudaReal sum(DeviceArray< cudaReal > const &in)
Compute sum of array elements (GPU kernel wrapper).
void eqV(DeviceArray< cudaReal > &a, DeviceArray< cudaReal > const &b, const int beginIdA, const int beginIdB, const int n)
Vector assignment, a[i] = b[i], kernel wrapper (cudaReal).
void sqSqNormV(DeviceArray< cudaReal > &a, DeviceArray< cudaComplex > const &b)
Norm of complex number to the 4th power, a[i] = norm(b[i])^4, kernel wrapper.
void mulEqV(DeviceArray< cudaReal > &a, DeviceArray< cudaReal > const &b, const int beginIdA, const int beginIdB, const int n)
Vector multiplication in-place, a[i] *= b[i], kernel wrapper (cudaReal).
void eqS(DeviceArray< cudaReal > &a, cudaReal const b, const int beginIdA, const int n)
Vector assignment, a[i] = b, kernel wrapper (cudaReal).
Periodic fields and crystallography.
PSCF package top-level namespace.
Utility classes for scientific computation.