Nektar++
ProcessInterpField.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: ProcessInterpField.cpp
4//
5// For more information, please see: http://www.nektar.info/
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Interpolate one field to another.
32//
33////////////////////////////////////////////////////////////////////////////////
34#include <iostream>
35#include <string>
36
37#include <boost/geometry.hpp>
38
43
44#include "ProcessInterpField.h"
45
46using namespace std;
47namespace bg = boost::geometry;
48namespace bgi = boost::geometry::index;
49
50namespace Nektar::FieldUtils
51{
52
56 "Interpolates one field to another, requires fromxml, "
57 "fromfld to be defined");
58
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}
76
78{
79}
80
81void ProcessInterpField::v_Process(po::variables_map &vm)
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();
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
250 if (m_f->m_verbose && m_f->m_comm->TreatAsRankZero())
251 {
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}
281
283 const int goal) const
284{
285 LibUtilities::PrintProgressbar(position, goal, "Interpolating");
286}
287} // namespace Nektar::FieldUtils
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
A class that contains algorithms for interpolation between pts fields, expansions and different meshe...
FIELD_UTILS_EXPORT void Interpolate(const T expInField, T &expOutField, NekDouble def_value=0.0)
Interpolate from an expansion to an expansion.
FieldSharedPtr m_f
Field object.
Definition: Module.h:239
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:272
void PrintProgressbar(const int position, const int goal) const
void v_Process(po::variables_map &vm) override
Write mesh to output file.
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
Abstract base class for processing modules.
Definition: Module.h:301
void SetProgressCallback(FuncPointerT func, ObjectPointerT obj)
sets a callback funtion which gets called every time the interpolation progresses
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static SessionReaderSharedPtr CreateInstance(int argc, char *argv[])
Creates an instance of the SessionReader class.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
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
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:47
int PrintProgressbar(const int position, const int goal, const std::string message, int lastprogress=-1)
Prints a progressbar.
Definition: Progressbar.hpp:65
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
Represents a command-line configuration option.
Definition: Module.h:129