Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::FieldUtils::ProcessInterpField Class Reference

This processing module interpolates one field to another. More...

#include <ProcessInterpField.h>

Inheritance diagram for Nektar::FieldUtils::ProcessInterpField:
[legend]

Public Member Functions

 ProcessInterpField (FieldSharedPtr f)
 
 ~ProcessInterpField () override
 
- Public Member Functions inherited from Nektar::FieldUtils::ProcessModule
 ProcessModule ()
 
 ProcessModule (FieldSharedPtr p_f)
 
- Public Member Functions inherited from Nektar::FieldUtils::Module
FIELD_UTILS_EXPORT Module (FieldSharedPtr p_f)
 
virtual ~Module ()=default
 
void Process (po::variables_map &vm)
 
std::string GetModuleName ()
 
std::string GetModuleDescription ()
 
const ConfigOptionGetConfigOption (const std::string &key) const
 
ModulePriority GetModulePriority ()
 
std::vector< ModuleKeyGetModulePrerequisites ()
 
FIELD_UTILS_EXPORT void RegisterConfig (std::string key, std::string value="")
 Register a configuration option with a module. More...
 
FIELD_UTILS_EXPORT void PrintConfig ()
 Print out all configuration options for a module. More...
 
FIELD_UTILS_EXPORT void SetDefaults ()
 Sets default configuration options for those which have not been set. More...
 
FIELD_UTILS_EXPORT void AddFile (std::string fileType, std::string fileName)
 
FIELD_UTILS_EXPORT void EvaluateTriFieldAtEquiSpacedPts (LocalRegions::ExpansionSharedPtr &exp, const Array< OneD, const NekDouble > &infield, Array< OneD, NekDouble > &outfield)
 

Static Public Member Functions

static std::shared_ptr< Modulecreate (FieldSharedPtr f)
 Creates an instance of this class. More...
 

Static Public Attributes

static ModuleKey className
 

Protected Member Functions

void v_Process (po::variables_map &vm) override
 Write mesh to output file. More...
 
void PrintProgressbar (const int position, const int goal) const
 
std::string v_GetModuleName () override
 
std::string v_GetModuleDescription () override
 
ModulePriority v_GetModulePriority () override
 
- Protected Member Functions inherited from Nektar::FieldUtils::Module
 Module ()
 
virtual void v_Process (po::variables_map &vm)
 
virtual std::string v_GetModuleName ()
 
virtual std::string v_GetModuleDescription ()
 
virtual ModulePriority v_GetModulePriority ()
 
virtual std::vector< ModuleKeyv_GetModulePrerequisites ()
 

Additional Inherited Members

- Public Attributes inherited from Nektar::FieldUtils::Module
FieldSharedPtr m_f
 Field object. More...
 
- Protected Attributes inherited from Nektar::FieldUtils::Module
std::map< std::string, ConfigOptionm_config
 List of configuration values. More...
 
std::set< std::string > m_allowedFiles
 List of allowed file formats. More...
 

Detailed Description

This processing module interpolates one field to another.

Definition at line 46 of file ProcessInterpField.h.

Constructor & Destructor Documentation

◆ ProcessInterpField()

Nektar::FieldUtils::ProcessInterpField::ProcessInterpField ( FieldSharedPtr  f)

Definition at line 59 of file ProcessInterpField.cpp.

59 : ProcessModule(f)
60{
61
62 m_config["fromxml"] = ConfigOption(
63 false, "NotSet", "Xml file from which to interpolate field");
64 m_config["fromfld"] = ConfigOption(
65 false, "NotSet", "Fld file from which to interpolate field");
66
67 m_config["clamptolowervalue"] =
68 ConfigOption(false, "-10000000", "Lower bound for interpolation value");
69 m_config["clamptouppervalue"] =
70 ConfigOption(false, "10000000", "Upper bound for interpolation value");
71 m_config["defaultvalue"] =
72 ConfigOption(false, "0", "Default value if point is outside domain");
73 m_config["realmodetoimag"] =
74 ConfigOption(false, "NotSet", "Take fields as sin mode");
75}
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:272

References Nektar::FieldUtils::Module::m_config.

◆ ~ProcessInterpField()

Nektar::FieldUtils::ProcessInterpField::~ProcessInterpField ( )
override

Definition at line 77 of file ProcessInterpField.cpp.

78{
79}

Member Function Documentation

◆ create()

static std::shared_ptr< Module > Nektar::FieldUtils::ProcessInterpField::create ( FieldSharedPtr  f)
inlinestatic

Creates an instance of this class.

Definition at line 50 of file ProcessInterpField.h.

51 {
53 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ PrintProgressbar()

void Nektar::FieldUtils::ProcessInterpField::PrintProgressbar ( const int  position,
const int  goal 
) const
protected

Definition at line 282 of file ProcessInterpField.cpp.

284{
285 LibUtilities::PrintProgressbar(position, goal, "Interpolating");
286}
int PrintProgressbar(const int position, const int goal, const std::string message, int lastprogress=-1)
Prints a progressbar.
Definition: Progressbar.hpp:65

References Nektar::LibUtilities::PrintProgressbar().

Referenced by v_Process().

◆ v_GetModuleDescription()

std::string Nektar::FieldUtils::ProcessInterpField::v_GetModuleDescription ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 70 of file ProcessInterpField.h.

71 {
72 return "Interpolating field";
73 }

◆ v_GetModuleName()

std::string Nektar::FieldUtils::ProcessInterpField::v_GetModuleName ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 65 of file ProcessInterpField.h.

66 {
67 return "ProcessInterpField";
68 }

◆ v_GetModulePriority()

ModulePriority Nektar::FieldUtils::ProcessInterpField::v_GetModulePriority ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 75 of file ProcessInterpField.h.

76 {
77 return eFillExp;
78 }

References Nektar::FieldUtils::eFillExp.

◆ v_Process()

void Nektar::FieldUtils::ProcessInterpField::v_Process ( po::variables_map &  vm)
overrideprotectedvirtual

Write mesh to output file.

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 81 of file ProcessInterpField.cpp.

82{
83 m_f->SetUpExp(vm);
84
85 FieldSharedPtr fromField = std::shared_ptr<Field>(new Field());
86
87 std::vector<std::string> files;
88
89 // set up session file for from field
90 char *argv[] = {const_cast<char *>("FieldConvert"), nullptr};
91 ParseUtils::GenerateVector(m_config["fromxml"].as<string>(), files);
93 1, argv, files,
94 LibUtilities::GetCommFactory().CreateInstance("Serial", 0, 0));
95
96 // Set up range based on min and max of local parallel partition
99
100 int numHomoDir = m_f->m_numHomogeneousDir;
101 int coordim = m_f->m_exp[0]->GetCoordim(0) + numHomoDir;
102 int npts = m_f->m_exp[0]->GetTotPoints();
103 Array<OneD, Array<OneD, NekDouble>> coords(3);
104
105 for (int i = 0; i < coordim; ++i)
106 {
107 coords[i] = Array<OneD, NekDouble>(npts);
108 }
109
110 for (int i = coordim; i < 3; ++i)
111 {
112 coords[i] = NullNekDouble1DArray;
113 }
114
115 m_f->m_exp[0]->GetCoords(coords[0], coords[1], coords[2]);
116
117 rng->m_checkShape = false;
118 switch (coordim)
119 {
120 case 3:
121 rng->m_doZrange = true;
122 rng->m_zmin = Vmath::Vmin(npts, coords[2], 1);
123 rng->m_zmax = Vmath::Vmax(npts, coords[2], 1);
124 /* Falls through. */
125 case 2:
126 rng->m_doYrange = true;
127 rng->m_ymin = Vmath::Vmin(npts, coords[1], 1);
128 rng->m_ymax = Vmath::Vmax(npts, coords[1], 1);
129 /* Falls through. */
130 case 1:
131 rng->m_doXrange = true;
132 rng->m_xmin = Vmath::Vmin(npts, coords[0], 1);
133 rng->m_xmax = Vmath::Vmax(npts, coords[0], 1);
134 break;
135 default:
136 NEKERROR(ErrorUtil::efatal, "coordim should be <= 3");
137 }
138
139 // setup rng parameters.
140 fromField->m_graph =
141 SpatialDomains::MeshGraph::Read(fromField->m_session, rng);
142
143 // Read in local from field partitions
144 const SpatialDomains::ExpansionInfoMap &expansions =
145 fromField->m_graph->GetExpansionInfo();
146
147 // check for case where no elements are specified on this
148 // parallel partition
149 if (!expansions.size())
150 {
151 return;
152 }
153
154 Array<OneD, int> ElementGIDs(expansions.size());
155
156 int i = 0;
157 for (auto &expIt : expansions)
158 {
159 ElementGIDs[i++] = expIt.second->m_geomShPtr->GetGlobalID();
160 }
161
162 string fromfld = m_config["fromfld"].as<string>();
163 m_f->FieldIOForFile(fromfld)->Import(
164 fromfld, fromField->m_fielddef, fromField->m_data,
166
167 int fromNumHomoDir = fromField->m_fielddef[0]->m_numHomogeneousDir;
168 for (i = 0; i < fromField->m_fielddef.size(); ++i)
169 {
170 int d1 = fromField->m_fielddef[i]->m_basis.size();
171 d1 -= 1;
172 if (d1 >= 0 && (fromField->m_fielddef[i]->m_basis[d1] ==
174 fromField->m_fielddef[i]->m_basis[d1] ==
176 {
177 fromField->m_fielddef[i]->m_homogeneousZIDs[0] += 2;
178 fromField->m_fielddef[i]->m_numModes[d1] = 4;
179 fromField->m_fielddef[i]->m_basis[d1] = LibUtilities::eFourier;
180 }
181 }
182
183 //----------------------------------------------
184 // Set up Expansion information to use mode order from field
185 fromField->m_graph->SetExpansionInfo(fromField->m_fielddef);
186
187 int nfields = fromField->m_fielddef[0]->m_fields.size();
188
189 fromField->m_exp.resize(nfields);
190 fromField->m_exp[0] = fromField->SetUpFirstExpList(fromNumHomoDir, true);
191
192 m_f->m_exp.resize(nfields);
193
194 // declare auxiliary fields.
195 for (i = 1; i < nfields; ++i)
196 {
197 m_f->m_exp[i] = m_f->AppendExpList(numHomoDir);
198 fromField->m_exp[i] = fromField->AppendExpList(fromNumHomoDir);
199 }
200
201 // load field into expansion in fromfield.
202 set<int> sinmode;
203 if (m_config["realmodetoimag"].as<string>().compare("NotSet"))
204 {
205 ParseUtils::GenerateVariableSet(m_config["realmodetoimag"].as<string>(),
206 m_f->m_variables, sinmode);
207 }
208 for (int j = 0; j < nfields; ++j)
209 {
210 for (i = 0; i < fromField->m_fielddef.size(); i++)
211 {
212 fromField->m_exp[j]->ExtractDataToCoeffs(
213 fromField->m_fielddef[i], fromField->m_data[i],
214 fromField->m_fielddef[0]->m_fields[j],
215 fromField->m_exp[j]->UpdateCoeffs());
216 }
217 if (fromNumHomoDir == 1)
218 {
219 fromField->m_exp[j]->SetWaveSpace(true);
220 if (sinmode.count(j))
221 {
222 int Ncoeff = fromField->m_exp[j]->GetPlane(2)->GetNcoeffs();
224 Ncoeff, -1., fromField->m_exp[j]->GetPlane(2)->GetCoeffs(),
225 1, fromField->m_exp[j]->GetPlane(3)->UpdateCoeffs(), 1);
226 Vmath::Zero(Ncoeff,
227 fromField->m_exp[j]->GetPlane(2)->UpdateCoeffs(),
228 1);
229 }
230 }
231 fromField->m_exp[j]->BwdTrans(fromField->m_exp[j]->GetCoeffs(),
232 fromField->m_exp[j]->UpdatePhys());
233 }
234
235 int nq1 = m_f->m_exp[0]->GetTotPoints();
236
237 NekDouble clamp_low = m_config["clamptolowervalue"].as<NekDouble>();
238 NekDouble clamp_up = m_config["clamptouppervalue"].as<NekDouble>();
239 NekDouble def_value = m_config["defaultvalue"].as<NekDouble>();
240
241 for (int i = 0; i < nfields; i++)
242 {
243 for (int j = 0; j < nq1; ++j)
244 {
245 m_f->m_exp[i]->UpdatePhys()[j] = def_value;
246 }
247 }
248
249 Interpolator<std::vector<MultiRegions::ExpListSharedPtr>> interp;
250 if (m_f->m_verbose && m_f->m_comm->TreatAsRankZero())
251 {
252 interp.SetProgressCallback(&ProcessInterpField::PrintProgressbar, this);
253 }
254
255 interp.Interpolate(fromField->m_exp, m_f->m_exp);
256
257 if (m_f->m_verbose && m_f->m_comm->TreatAsRankZero())
258 {
259 cout << endl;
260 }
261
262 for (int i = 0; i < nfields; ++i)
263 {
264 for (int j = 0; j < nq1; ++j)
265 {
266 if (m_f->m_exp[i]->GetPhys()[j] > clamp_up)
267 {
268 m_f->m_exp[i]->UpdatePhys()[j] = clamp_up;
269 }
270 else if (m_f->m_exp[i]->GetPhys()[j] < clamp_low)
271 {
272 m_f->m_exp[i]->UpdatePhys()[j] = clamp_low;
273 }
274 }
275 m_f->m_exp[i]->FwdTransLocalElmt(m_f->m_exp[i]->GetPhys(),
276 m_f->m_exp[i]->UpdateCoeffs());
277 }
278 // save field names
279 m_f->m_variables = fromField->m_fielddef[0]->m_fields;
280}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
FieldSharedPtr m_f
Field object.
Definition: Module.h:239
void PrintProgressbar(const int position, const int goal) const
static SessionReaderSharedPtr CreateInstance(int argc, char *argv[])
Creates an instance of the SessionReader class.
static bool GenerateVariableSet(const std::string &str, const std::vector< std::string > &variables, std::set< int > &out)
Generate a set of variable locations.
Definition: ParseUtils.cpp:168
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:130
static MeshGraphSharedPtr Read(const LibUtilities::SessionReaderSharedPtr pSession, LibUtilities::DomainRangeShPtr rng=LibUtilities::NullDomainRangeShPtr, bool fillGraph=true, SpatialDomains::MeshGraphSharedPtr partitionedGraph=nullptr)
Definition: MeshGraph.cpp:115
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:990
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:51
std::shared_ptr< DomainRange > DomainRangeShPtr
Definition: DomainRange.h:64
CommFactory & GetCommFactory()
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode .
Definition: BasisType.h:68
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode .
Definition: BasisType.h:66
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
Definition: MeshGraph.h:141
static Array< OneD, NekDouble > NullNekDouble1DArray
double NekDouble
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
Definition: Vmath.hpp:725
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
Definition: Vmath.hpp:644

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::SessionReader::CreateInstance(), Nektar::ErrorUtil::efatal, Nektar::LibUtilities::eFourier, Nektar::LibUtilities::eFourierHalfModeIm, Nektar::LibUtilities::eFourierHalfModeRe, Nektar::ParseUtils::GenerateVariableSet(), Nektar::ParseUtils::GenerateVector(), Nektar::LibUtilities::GetCommFactory(), Nektar::FieldUtils::Interpolator< T >::Interpolate(), Nektar::FieldUtils::Module::m_config, Nektar::FieldUtils::Module::m_f, NEKERROR, Nektar::LibUtilities::NullFieldMetaDataMap, Nektar::NullNekDouble1DArray, PrintProgressbar(), Nektar::SpatialDomains::MeshGraph::Read(), Nektar::LibUtilities::Interpolator::SetProgressCallback(), Vmath::Smul(), Vmath::Vmax(), Vmath::Vmin(), and Vmath::Zero().

Member Data Documentation

◆ className

ModuleKey Nektar::FieldUtils::ProcessInterpField::className
static
Initial value:
=
"Interpolates one field to another, requires fromxml, "
"fromfld to be defined")
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:47

Definition at line 54 of file ProcessInterpField.h.