Nektar++
Public Member Functions | Public Attributes | Private Attributes | List of all members
Nektar::FieldUtils::Field Struct Reference

#include <Field.hpp>

Public Member Functions

FIELD_UTILS_EXPORT Field ()
 
FIELD_UTILS_EXPORT ~Field ()
 
FIELD_UTILS_EXPORT void SetUpExp (boost::program_options::variables_map &vm)
 
FIELD_UTILS_EXPORT void CreateExp (boost::program_options::variables_map &vm, bool newExp)
 
FIELD_UTILS_EXPORT MultiRegions::ExpListSharedPtr SetUpFirstExpList (int NumHomogeneousDir, bool fldfilegiven=false)
 
FIELD_UTILS_EXPORT LibUtilities::FieldIOSharedPtr FieldIOForFile (std::string filename)
 Construct a FieldIO object for the file filename. More...
 
FIELD_UTILS_EXPORT MultiRegions::ExpListSharedPtr AppendExpList (int NumHomogeneousDir, std::string var="DefaultVar", bool NewField=false)
 
FIELD_UTILS_EXPORT void ClearField ()
 

Public Attributes

bool m_verbose
 
std::vector< LibUtilities::FieldDefinitionsSharedPtrm_fielddef
 
std::vector< std::vector< double > > m_data
 
std::vector< MultiRegions::ExpListSharedPtrm_exp
 
std::vector< std::string > m_variables
 
int m_numHomogeneousDir
 
bool m_declareExpansionAsContField
 
bool m_declareExpansionAsDisContField
 
bool m_requireBoundaryExpansion
 
bool m_useFFT
 
LibUtilities::CommSharedPtr m_comm
 
LibUtilities::CommSharedPtr m_defComm
 
LibUtilities::CommSharedPtr m_partComm
 
int m_nParts = 1
 
po::variables_map m_vm
 
LibUtilities::SessionReaderSharedPtr m_session
 
SpatialDomains::MeshGraphSharedPtr m_graph
 
std::map< std::string, std::vector< std::string > > m_inputfiles
 
bool m_writeBndFld
 
std::vector< unsigned int > m_bndRegionsToWrite
 
bool m_addNormals
 
LibUtilities::PtsFieldSharedPtr m_fieldPts
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 

Private Attributes

std::map< std::string, LibUtilities::FieldIOSharedPtrm_fld
 Map to store FieldIO instances. Key is the reader type, value is the FieldIO object. More...
 

Detailed Description

Definition at line 62 of file Field.hpp.

Constructor & Destructor Documentation

◆ Field()

FIELD_UTILS_EXPORT Nektar::FieldUtils::Field::Field ( )
inline

Definition at line 64 of file Field.hpp.

69  {
70  }
static PtsFieldSharedPtr NullPtsField
Definition: PtsField.h:191
bool m_declareExpansionAsDisContField
Definition: Field.hpp:88
bool m_declareExpansionAsContField
Definition: Field.hpp:87
LibUtilities::PtsFieldSharedPtr m_fieldPts
Definition: Field.hpp:108

◆ ~Field()

FIELD_UTILS_EXPORT Nektar::FieldUtils::Field::~Field ( )
inline

Definition at line 72 of file Field.hpp.

73  {
74  if (m_comm)
75  {
76  m_comm->Finalise();
77  }
78  }
LibUtilities::CommSharedPtr m_comm
Definition: Field.hpp:94

References m_comm.

Member Function Documentation

◆ AppendExpList()

FIELD_UTILS_EXPORT MultiRegions::ExpListSharedPtr Nektar::FieldUtils::Field::AppendExpList ( int  NumHomogeneousDir,
std::string  var = "DefaultVar",
bool  NewField = false 
)
inline

Definition at line 662 of file Field.hpp.

665  {
666  if (var.compare("DefaultVar") == 0 && m_requireBoundaryExpansion)
667  {
668  if (m_session->GetVariables().size())
669  {
670  var = m_session->GetVariables()[0];
671  }
672  }
674  switch (m_graph->GetMeshDimension())
675  {
676  case 1:
677  {
678  if (NumHomogeneousDir == 1)
679  {
682  "ContFieldHomogeneous1D or "
683  "DisContFieldHomogenenous1D has not been "
684  "implemented");
685 
687  std::dynamic_pointer_cast<
688  MultiRegions::ExpList2DHomogeneous1D>(m_exp[0]);
689 
691  AllocateSharedPtr(*tmp2);
692  }
693  else if (NumHomogeneousDir == 2)
694  {
696  {
698  std::dynamic_pointer_cast<
699  MultiRegions::ContField3DHomogeneous2D>(
700  m_exp[0]);
701 
702  tmp = MemoryManager<
703  MultiRegions::ContField3DHomogeneous2D>::
704  AllocateSharedPtr(*tmp2);
705  }
707  {
709  tmp2 = std::dynamic_pointer_cast<
710  MultiRegions::DisContField3DHomogeneous2D>(
711  m_exp[0]);
712 
713  tmp = MemoryManager<
714  MultiRegions::DisContField3DHomogeneous2D>::
715  AllocateSharedPtr(*tmp2);
716  }
717  else
718  {
720  std::dynamic_pointer_cast<
721  MultiRegions::ExpList3DHomogeneous2D>(m_exp[0]);
722 
723  tmp = MemoryManager<
724  MultiRegions::ExpList3DHomogeneous2D>::
725  AllocateSharedPtr(*tmp2);
726  }
727  }
728  else
729  {
731  {
733  std::dynamic_pointer_cast<MultiRegions::ContField>(
734  m_exp[0]);
735 
738  false,
740  }
742  {
744  std::dynamic_pointer_cast<
745  MultiRegions::DisContField>(m_exp[0]);
746 
750  }
751  else
752  {
754  std::dynamic_pointer_cast<MultiRegions::ExpList>(
755  m_exp[0]);
756 
757  tmp = MemoryManager<
758  MultiRegions::ExpList>::AllocateSharedPtr(*tmp2);
759  }
760  }
761  }
762  break;
763  case 2:
764  {
765  if (NumHomogeneousDir == 1)
766  {
768  {
769  if (NewField)
770  {
771  bool dealiasing = false;
772 
773  tmp = MemoryManager<
774  MultiRegions::ContField3DHomogeneous1D>::
775  AllocateSharedPtr(m_session,
776  m_exp[0]
777  ->GetHomogeneousBasis()
778  ->GetBasisKey(),
779  m_exp[0]->GetHomoLen(),
780  m_useFFT, dealiasing, m_graph,
781  var,
783  }
784  else
785  {
787  tmp2 = std::dynamic_pointer_cast<
788  MultiRegions::ContField3DHomogeneous1D>(
789  m_exp[0]);
790 
791  ASSERTL0(tmp2, "Failed to type cast m_exp[0]");
792  tmp = MemoryManager<
793  MultiRegions::ContField3DHomogeneous1D>::
794  AllocateSharedPtr(*tmp2, m_graph, var);
795  }
796  }
798  {
799  if (NewField)
800  {
801  bool dealiasing = false;
802 
803  tmp = MemoryManager<
804  MultiRegions::DisContField3DHomogeneous1D>::
805  AllocateSharedPtr(m_session,
806  m_exp[0]
807  ->GetHomogeneousBasis()
808  ->GetBasisKey(),
809  m_exp[0]->GetHomoLen(),
810  m_useFFT, dealiasing, m_graph,
811  var,
813  }
814  else
815  {
817  tmp2 = std::dynamic_pointer_cast<
818  MultiRegions::DisContField3DHomogeneous1D>(
819  m_exp[0]);
820  ASSERTL0(tmp2, "Failed to type cast m_exp[0]");
821 
822  tmp = MemoryManager<
823  MultiRegions::DisContField3DHomogeneous1D>::
824  AllocateSharedPtr(*tmp2);
825  }
826  }
827  else
828  {
829  if (NewField)
830  {
831  bool dealiasing = false;
832 
833  tmp = MemoryManager<
834  MultiRegions::ExpList3DHomogeneous1D>::
835  AllocateSharedPtr(m_session,
836  m_exp[0]
837  ->GetHomogeneousBasis()
838  ->GetBasisKey(),
839  m_exp[0]->GetHomoLen(),
840  m_useFFT, dealiasing, m_graph,
841  "DefaultVar",
843  }
844  else
845  {
847  std::dynamic_pointer_cast<
848  MultiRegions::ExpList3DHomogeneous1D>(
849  m_exp[0]);
850  ASSERTL0(tmp2, "Failed to type cast m_exp[0]");
851 
852  tmp = MemoryManager<
853  MultiRegions::ExpList3DHomogeneous1D>::
854  AllocateSharedPtr(*tmp2);
855  }
856  }
857  }
858  else
859  {
861  {
862  if (NewField)
863  {
866  false,
868  }
869  else // call copy constructor
870  {
871 
873  std::dynamic_pointer_cast<
874  MultiRegions::ContField>(m_exp[0]);
875 
877  AllocateSharedPtr(*tmp2, m_graph, var);
878  }
879  }
881  {
882  if (NewField)
883  {
886  true,
888  }
889  else // call copy constructor
890  {
892  std::dynamic_pointer_cast<
893  MultiRegions::DisContField>(m_exp[0]);
894 
896  AllocateSharedPtr(*tmp2, m_graph, var);
897  }
898  }
899  else
900  {
902  std::dynamic_pointer_cast<MultiRegions::ExpList>(
903  m_exp[0]);
904 
905  tmp = MemoryManager<
906  MultiRegions::ExpList>::AllocateSharedPtr(*tmp2);
907  }
908  }
909  }
910  break;
911  case 3:
912  {
914  {
915  if (NewField)
916  {
919  false,
921  }
922  else
923  {
925  std::dynamic_pointer_cast<MultiRegions::ContField>(
926  m_exp[0]);
927 
928  tmp = MemoryManager<
929  MultiRegions::ContField>::AllocateSharedPtr(*tmp2,
930  m_graph,
931  var);
932  }
933  }
935  {
936  if (NewField)
937  {
941  }
942  else
943  {
945  std::dynamic_pointer_cast<
946  MultiRegions::DisContField>(m_exp[0]);
947 
949  AllocateSharedPtr(*tmp2, m_graph, var);
950  }
951  }
952  else
953  {
955  std::dynamic_pointer_cast<MultiRegions::ExpList>(
956  m_exp[0]);
957 
958  tmp =
960  *tmp2);
961  }
962  }
963  break;
964  default:
965  ASSERTL0(false, "Expansion dimension not recognised");
966  break;
967  }
968 
969  return tmp;
970  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< DisContField3DHomogeneous2D > DisContField3DHomogeneous2DSharedPtr
std::shared_ptr< DisContField3DHomogeneous1D > DisContField3DHomogeneous1DSharedPtr
std::shared_ptr< DisContField > DisContFieldSharedPtr
Definition: DisContField.h:377
std::shared_ptr< ContField3DHomogeneous2D > ContField3DHomogeneous2DSharedPtr
std::shared_ptr< ContField3DHomogeneous1D > ContField3DHomogeneous1DSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< ExpList2DHomogeneous1D > ExpList2DHomogeneous1DSharedPtr
Shared pointer to an ExpList2DHomogeneous1D object.
std::shared_ptr< ExpList3DHomogeneous2D > ExpList3DHomogeneous2DSharedPtr
Shared pointer to an ExpList3DHomogeneous2D object.
std::shared_ptr< ExpList3DHomogeneous1D > ExpList3DHomogeneous1DSharedPtr
Shared pointer to an ExpList3DHomogeneous1D object.
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:277
SpatialDomains::MeshGraphSharedPtr m_graph
Definition: Field.hpp:101
LibUtilities::SessionReaderSharedPtr m_session
Definition: Field.hpp:100
std::vector< MultiRegions::ExpListSharedPtr > m_exp
Definition: Field.hpp:82

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::Collections::eNoCollection, m_declareExpansionAsContField, m_declareExpansionAsDisContField, m_exp, m_graph, m_requireBoundaryExpansion, m_session, and m_useFFT.

Referenced by CreateExp().

◆ ClearField()

FIELD_UTILS_EXPORT void Nektar::FieldUtils::Field::ClearField ( )
inline

Definition at line 972 of file Field.hpp.

973  {
977  m_exp.clear();
978  m_fielddef = std::vector<LibUtilities::FieldDefinitionsSharedPtr>();
979  m_data = std::vector<std::vector<NekDouble>>();
980  m_variables.clear();
981  }
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
std::vector< std::string > m_variables
Definition: Field.hpp:83
std::vector< LibUtilities::FieldDefinitionsSharedPtr > m_fielddef
Definition: Field.hpp:80
std::vector< std::vector< double > > m_data
Definition: Field.hpp:81

References m_data, m_exp, m_fielddef, m_fieldPts, m_graph, m_session, m_variables, and Nektar::LibUtilities::NullPtsField.

Referenced by export_Field().

◆ CreateExp()

FIELD_UTILS_EXPORT void Nektar::FieldUtils::Field::CreateExp ( boost::program_options::variables_map &  vm,
bool  newExp 
)
inline

Definition at line 124 of file Field.hpp.

126  {
127 
128  bool fldfilegiven = (m_fielddef.size() != 0);
129  if (newExp)
130  {
131  LibUtilities::Timer timerpart;
132  if (m_verbose)
133  {
134  if (m_comm->TreatAsRankZero())
135  {
136  timerpart.Start();
137  }
138  }
139  // check to see if fld file defined so can use in
140  // expansion defintion if required
141 
142  bool expFromFld = fldfilegiven && !vm.count("useSessionExpansion");
143 
144  // load fielddef header if fld file is defined. This gives
145  // precedence to Homogeneous definition in fld file
147  if (expFromFld)
148  {
149  m_numHomogeneousDir = m_fielddef[0]->m_numHomogeneousDir;
150 
151  // Set up Expansion information to use mode order from field
152  m_graph->SetExpansionInfo(m_fielddef);
153  }
154  else
155  {
156  if (m_session->DefinesSolverInfo("HOMOGENEOUS"))
157  {
158  std::string HomoStr =
159  m_session->GetSolverInfo("HOMOGENEOUS");
160 
161  if ((HomoStr == "HOMOGENEOUS1D") ||
162  (HomoStr == "Homogeneous1D") || (HomoStr == "1D") ||
163  (HomoStr == "Homo1D"))
164  {
166  }
167  if ((HomoStr == "HOMOGENEOUS2D") ||
168  (HomoStr == "Homogeneous2D") || (HomoStr == "2D") ||
169  (HomoStr == "Homo2D"))
170  {
172  }
173  }
174  }
175 
176  m_exp.resize(1);
177  // Check if there are any elements to process
178  std::vector<int> IDs;
179  auto domain = m_graph->GetDomain();
180  for (int d = 0; d < domain.size(); ++d)
181  {
182  for (auto &compIter : domain[d])
183  {
184  for (auto &x : compIter.second->m_geomVec)
185  {
186  IDs.push_back(x->GetGlobalID());
187  }
188  }
189  }
190  // if Range has been specified it is possible to have a
191  // partition which is empty so check this and return with empty
192  // expansion if no elements present.
193  if (!IDs.size())
194  {
195  m_exp[0] =
197  return;
198  }
199 
200  // Adjust number of quadrature points
201  if (vm.count("output-points"))
202  {
203  int nPointsNew = vm["output-points"].as<int>();
204  m_graph->SetExpansionInfoToPointOrder(nPointsNew);
205  }
206 
207  if (m_verbose)
208  {
209  if (m_comm->TreatAsRankZero())
210  {
211  timerpart.Stop();
212  NekDouble cpuTime = timerpart.TimePerTest(1);
213 
214  std::stringstream ss;
215  ss << cpuTime << "s";
216  std::cout << "\t CreateExp setexpansion CPU Time: "
217  << std::setw(8) << std::left << ss.str()
218  << std::endl;
219  timerpart.Start();
220  }
221  }
222  // Override number of planes with value from cmd line
223  if (m_numHomogeneousDir == 1 && vm.count("output-points-hom-z"))
224  {
225  int expdim = m_graph->GetMeshDimension();
226  m_fielddef[0]->m_numModes[expdim] =
227  vm["output-points-hom-z"].as<int>();
228  }
229  m_exp[0] = SetUpFirstExpList(m_numHomogeneousDir, expFromFld);
230  if (m_verbose)
231  {
232  if (m_comm->TreatAsRankZero())
233  {
234  timerpart.Stop();
235  NekDouble cpuTime = timerpart.TimePerTest(1);
236 
237  std::stringstream ss1;
238 
239  ss1 << cpuTime << "s";
240  std::cout << "\t CreateExp set first exp CPU Time: "
241  << std::setw(8) << std::left << ss1.str()
242  << std::endl;
243  }
244  }
245  }
246 
247  if (fldfilegiven)
248  {
249  int i, j, nfields, nstrips;
250  m_session->LoadParameter("Strip_Z", nstrips, 1);
251  std::vector<std::string> vars = m_session->GetVariables();
252 
253  if (vm.count("useSessionVariables"))
254  {
255  m_variables = vars;
256  }
257  nfields = m_variables.size();
258 
259  m_exp.resize(nfields * nstrips);
260  // declare other fields;
261  for (int s = 0; s < nstrips; ++s) // homogeneous strip varient
262  {
263  for (i = 0; i < nfields; ++i)
264  {
265  if (i < vars.size())
266  {
267  // check to see if field already defined
268  if (!m_exp[s * nfields + i])
269  {
270  m_exp[s * nfields + i] =
272  }
273  }
274  else
275  {
276  if (vars.size())
277  {
278  m_exp[s * nfields + i] =
280  }
281  else
282  {
283  m_exp[s * nfields + i] =
285  }
286  }
287  }
288  }
289 
290  // Extract data to coeffs and bwd transform
291  for (int s = 0; s < nstrips; ++s) // homogeneous strip varient
292  {
293  for (j = 0; j < nfields; ++j)
294  {
295  for (i = 0; i < m_data.size() / nstrips; ++i)
296  {
297  int n = i * nstrips + s;
298  // In case of multiple flds, we might not have a
299  // variable in this m_data[n] -> skip in this case
300  auto it =
301  find(m_fielddef[n]->m_fields.begin(),
302  m_fielddef[n]->m_fields.end(), m_variables[j]);
303  if (it != m_fielddef[n]->m_fields.end())
304  {
305  m_exp[s * nfields + j]->ExtractDataToCoeffs(
306  m_fielddef[n], m_data[n], m_variables[j],
307  m_exp[s * nfields + j]->UpdateCoeffs());
308  }
309  }
310  m_exp[s * nfields + j]->BwdTrans(
311  m_exp[s * nfields + j]->GetCoeffs(),
312  m_exp[s * nfields + j]->UpdatePhys());
313  }
314  }
315  // Clear fielddef and data
316  // (they should not be used after running this module)
317  m_fielddef = std::vector<LibUtilities::FieldDefinitionsSharedPtr>();
318  m_data = std::vector<std::vector<NekDouble>>();
319  }
320  }
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:327
double NekDouble
FIELD_UTILS_EXPORT MultiRegions::ExpListSharedPtr SetUpFirstExpList(int NumHomogeneousDir, bool fldfilegiven=false)
Definition: Field.hpp:322
FIELD_UTILS_EXPORT MultiRegions::ExpListSharedPtr AppendExpList(int NumHomogeneousDir, std::string var="DefaultVar", bool NewField=false)
Definition: Field.hpp:662

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), AppendExpList(), Nektar::StdRegions::find(), m_comm, m_data, m_exp, m_fielddef, m_graph, m_numHomogeneousDir, m_session, m_variables, m_verbose, SetUpFirstExpList(), Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and Nektar::LibUtilities::Timer::TimePerTest().

Referenced by SetUpExp().

◆ FieldIOForFile()

FIELD_UTILS_EXPORT LibUtilities::FieldIOSharedPtr Nektar::FieldUtils::Field::FieldIOForFile ( std::string  filename)
inline

Construct a FieldIO object for the file filename.

This routine constructs an appropriate FieldIO object for a filename through the LibUtilities::FieldIO::GetFileType function to detect the file format. The result is then cached in Field::m_fld to avoid needing to repeatedly construct the object.

Parameters
filenameFilename to open.
Returns
Reader for filename.

Definition at line 641 of file Field.hpp.

643  {
645  std::string fmt = LibUtilities::FieldIO::GetFileType(filename, c);
646 
647  auto it = m_fld.find(fmt);
648 
649  if (it == m_fld.end())
650  {
653  m_fld[fmt] = fld;
654  return fld;
655  }
656  else
657  {
658  return it->second;
659  }
660  }
static const std::string GetFileType(const std::string &filename, CommSharedPtr comm)
Determine file type of given input file.
Definition: FieldIO.cpp:97
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:301
FieldIOFactory & GetFieldIOFactory()
Returns the FieldIO factory.
Definition: FieldIO.cpp:72
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:54
std::map< std::string, LibUtilities::FieldIOSharedPtr > m_fld
Map to store FieldIO instances. Key is the reader type, value is the FieldIO object.
Definition: Field.hpp:986

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::GetFieldIOFactory(), Nektar::LibUtilities::FieldIO::GetFileType(), m_comm, and m_fld.

◆ SetUpExp()

FIELD_UTILS_EXPORT void Nektar::FieldUtils::Field::SetUpExp ( boost::program_options::variables_map &  vm)
inline

Definition at line 112 of file Field.hpp.

113  {
114  if (m_graph && !m_exp.size())
115  {
116  CreateExp(vm, true);
117  }
118  else if (m_graph && m_data.size())
119  {
120  CreateExp(vm, false);
121  }
122  }
FIELD_UTILS_EXPORT void CreateExp(boost::program_options::variables_map &vm, bool newExp)
Definition: Field.hpp:124

References CreateExp(), m_data, m_exp, and m_graph.

◆ SetUpFirstExpList()

FIELD_UTILS_EXPORT MultiRegions::ExpListSharedPtr Nektar::FieldUtils::Field::SetUpFirstExpList ( int  NumHomogeneousDir,
bool  fldfilegiven = false 
)
inline

Definition at line 322 of file Field.hpp.

324  {
325 
327 
328  // Set up expansion list
329  int expdim = m_graph->GetMeshDimension();
330  bool dealiasing = false;
331 
332  m_session->MatchSolverInfo("USEFFT", "FFTW", m_useFFT, false);
333 
334  switch (expdim)
335  {
336  case 1:
337  {
338  ASSERTL0(NumHomogeneousDir <= 2,
339  "Quasi-3D approach is only set up for 1 or 2 "
340  "homogeneous directions");
341 
342  if (NumHomogeneousDir == 1)
343  {
345 
346  // Define Homogeneous expansion
347  int nplanes;
348  NekDouble ly;
350 
351  if (fldfilegiven)
352  {
353  nplanes = m_fielddef[0]->m_numModes[1];
354  ly = m_fielddef[0]->m_homogeneousLengths[0];
355  btype = m_fielddef[0]->m_basis[1];
356  }
357  else
358  {
359  m_session->LoadParameter("HomModesZ", nplanes);
360  m_session->LoadParameter("LY", ly);
361  btype = LibUtilities::eFourier;
362  }
363 
364  // Choose points to be at evenly spaced points at
365  // nplanes points
366  const LibUtilities::PointsKey Pkey(
368 
369  const LibUtilities::BasisKey Bkey(btype, nplanes, Pkey);
370 
373  {
374  ASSERTL0(false, "ContFieldHomogeneous1D or "
375  "DisContFieldHomogenenous1D has "
376  "not been implemented");
377  }
378 
379  Exp2DH1 =
382  dealiasing, m_graph,
384  exp = Exp2DH1;
385  }
386  else if (NumHomogeneousDir == 2)
387  {
389 
390  int nylines, nzlines;
391  NekDouble ly, lz;
392  LibUtilities::BasisType btype1, btype2;
393 
394  if (fldfilegiven)
395  {
396  nylines = m_fielddef[0]->m_numModes[1];
397  nzlines = m_fielddef[0]->m_numModes[2];
398  ly = m_fielddef[0]->m_homogeneousLengths[0];
399  lz = m_fielddef[0]->m_homogeneousLengths[1];
400  btype1 = m_fielddef[0]->m_basis[1];
401  btype2 = m_fielddef[0]->m_basis[2];
402  }
403  else
404  {
405  m_session->LoadParameter("HomModesY", nylines);
406  m_session->LoadParameter("HomModesZ", nzlines);
407  m_session->LoadParameter("LY", ly);
408  m_session->LoadParameter("LZ", lz);
409  btype1 = LibUtilities::eFourier;
410  btype2 = LibUtilities::eFourier;
411  }
412 
413  // Choose points to be at evenly spaced points at
414  // nplanes points
415  const LibUtilities::PointsKey PkeyY(
417  const LibUtilities::BasisKey BkeyY(btype1, nylines, PkeyY);
418 
419  const LibUtilities::PointsKey PkeyZ(
421  const LibUtilities::BasisKey BkeyZ(btype2, nzlines, PkeyZ);
422 
424  {
425  Exp3DH2 = MemoryManager<
426  MultiRegions::ContField3DHomogeneous2D>::
427  AllocateSharedPtr(m_session, BkeyY, BkeyZ, ly, lz,
428  m_useFFT, dealiasing, m_graph,
429  m_session->GetVariable(0),
431  }
433  {
434  Exp3DH2 = MemoryManager<
435  MultiRegions::DisContField3DHomogeneous2D>::
436  AllocateSharedPtr(m_session, BkeyY, BkeyZ, ly, lz,
437  m_useFFT, dealiasing, m_graph,
438  m_session->GetVariable(0),
440  }
441  else
442  {
443  Exp3DH2 = MemoryManager<
444  MultiRegions::ExpList3DHomogeneous2D>::
445  AllocateSharedPtr(m_session, BkeyY, BkeyZ, ly, lz,
446  m_useFFT, dealiasing, m_graph,
448  }
449 
450  exp = Exp3DH2;
451  }
452  else
453  {
455 
457  {
460  m_session, m_graph, m_session->GetVariable(0),
461  true, false, Collections::eNoCollection);
462  }
464  {
467  m_session->GetVariable(0), true,
469  }
470  else
471  {
474  "DefaultVar",
476  }
477 
478  exp = Exp1D;
479  }
480  }
481  break;
482  case 2:
483  {
484  ASSERTL0(NumHomogeneousDir <= 1,
485  "NumHomogeneousDir is only set up for 1");
486 
487  if (NumHomogeneousDir == 1)
488  {
490 
491  // Define Homogeneous expansion
492  int nplanes;
493  NekDouble lz;
497 
498  if (fldfilegiven)
499  {
500  nplanes = m_fielddef[0]->m_numModes[2];
501  lz = m_fielddef[0]->m_homogeneousLengths[0];
502  btype = m_fielddef[0]->m_basis[2];
503 
505  {
506  btype = LibUtilities::eFourier;
507  m_fielddef[0]->m_basis[2] =
509  if (nplanes <= 2)
510  {
511  nplanes = 4;
512  }
513  }
514  else if (btype == LibUtilities::eFourierHalfModeRe &&
515  nplanes == 1)
516  {
518  }
519  }
520  else
521  {
522  m_session->LoadParameter("HomModesZ", nplanes);
523  m_session->LoadParameter("LZ", lz);
524  btype = LibUtilities::eFourier;
525  }
526  // Choose points to be at evenly spaced points at
527  // nplanes points
528  const LibUtilities::PointsKey Pkey(nplanes, ptype);
529 
530  const LibUtilities::BasisKey Bkey(btype, nplanes, Pkey);
531 
533  {
534  Exp3DH1 = MemoryManager<
535  MultiRegions::ContField3DHomogeneous1D>::
536  AllocateSharedPtr(m_session, Bkey, lz, m_useFFT,
537  dealiasing, m_graph,
538  m_session->GetVariable(0),
540  }
541 
543  {
544  Exp3DH1 = MemoryManager<
545  MultiRegions::DisContField3DHomogeneous1D>::
546  AllocateSharedPtr(m_session, Bkey, lz, m_useFFT,
547  dealiasing, m_graph,
548  m_session->GetVariable(0),
550  }
551 
552  else
553  {
554  Exp3DH1 = MemoryManager<
555  MultiRegions::ExpList3DHomogeneous1D>::
556  AllocateSharedPtr(m_session, Bkey, lz, m_useFFT,
557  dealiasing, m_graph, "DefaultVar",
559  }
560  exp = Exp3DH1;
561  }
562  else
563  {
565 
567  {
570  m_session, m_graph, m_session->GetVariable(0),
571  true, false, Collections::eNoCollection);
572  }
573 
575  {
578  m_session->GetVariable(0), true,
580  }
581  else
582  {
585  "DefaultVar",
587  }
588 
589  exp = Exp2D;
590  }
591  }
592  break;
593  case 3:
594  {
596 
598  {
601  m_session->GetVariable(0), true,
603  }
605  {
608  m_session->GetVariable(0), true,
610  }
611  else
612  {
613  Exp3D =
615  m_session, m_graph, true, "DefaultVar",
617  }
618 
619  exp = Exp3D;
620  }
621  break;
622  default:
623  ASSERTL0(false, "Expansion dimension not recognised");
624  break;
625  }
626 
627  return exp;
628  };
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:76
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:75
@ eFourierSingleMode
Fourier ModifiedExpansion with just the first mode .
Definition: BasisType.h:66
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode .
Definition: BasisType.h:68
@ eFourier
Fourier Expansion .
Definition: BasisType.h:57

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::LibUtilities::eFourier, Nektar::LibUtilities::eFourierEvenlySpaced, Nektar::LibUtilities::eFourierHalfModeRe, Nektar::LibUtilities::eFourierSingleMode, Nektar::Collections::eNoCollection, Nektar::LibUtilities::ePolyEvenlySpaced, m_declareExpansionAsContField, m_declareExpansionAsDisContField, m_fielddef, m_graph, m_session, and m_useFFT.

Referenced by CreateExp().

Member Data Documentation

◆ m_addNormals

bool Nektar::FieldUtils::Field::m_addNormals

Definition at line 106 of file Field.hpp.

◆ m_bndRegionsToWrite

std::vector<unsigned int> Nektar::FieldUtils::Field::m_bndRegionsToWrite

Definition at line 105 of file Field.hpp.

◆ m_comm

LibUtilities::CommSharedPtr Nektar::FieldUtils::Field::m_comm

Definition at line 94 of file Field.hpp.

Referenced by CreateExp(), FieldIOForFile(), and ~Field().

◆ m_data

std::vector<std::vector<double> > Nektar::FieldUtils::Field::m_data

Definition at line 81 of file Field.hpp.

Referenced by ClearField(), CreateExp(), and SetUpExp().

◆ m_declareExpansionAsContField

bool Nektar::FieldUtils::Field::m_declareExpansionAsContField

Definition at line 87 of file Field.hpp.

Referenced by AppendExpList(), and SetUpFirstExpList().

◆ m_declareExpansionAsDisContField

bool Nektar::FieldUtils::Field::m_declareExpansionAsDisContField

Definition at line 88 of file Field.hpp.

Referenced by AppendExpList(), and SetUpFirstExpList().

◆ m_defComm

LibUtilities::CommSharedPtr Nektar::FieldUtils::Field::m_defComm

Definition at line 95 of file Field.hpp.

◆ m_exp

std::vector<MultiRegions::ExpListSharedPtr> Nektar::FieldUtils::Field::m_exp

Definition at line 82 of file Field.hpp.

Referenced by AppendExpList(), ClearField(), CreateExp(), and SetUpExp().

◆ m_fielddef

std::vector<LibUtilities::FieldDefinitionsSharedPtr> Nektar::FieldUtils::Field::m_fielddef

Definition at line 80 of file Field.hpp.

Referenced by ClearField(), CreateExp(), and SetUpFirstExpList().

◆ m_fieldMetaDataMap

LibUtilities::FieldMetaDataMap Nektar::FieldUtils::Field::m_fieldMetaDataMap

Definition at line 110 of file Field.hpp.

◆ m_fieldPts

LibUtilities::PtsFieldSharedPtr Nektar::FieldUtils::Field::m_fieldPts

Definition at line 108 of file Field.hpp.

Referenced by ClearField().

◆ m_fld

std::map<std::string, LibUtilities::FieldIOSharedPtr> Nektar::FieldUtils::Field::m_fld
private

Map to store FieldIO instances. Key is the reader type, value is the FieldIO object.

Definition at line 986 of file Field.hpp.

Referenced by FieldIOForFile().

◆ m_graph

SpatialDomains::MeshGraphSharedPtr Nektar::FieldUtils::Field::m_graph

Definition at line 101 of file Field.hpp.

Referenced by AppendExpList(), ClearField(), CreateExp(), SetUpExp(), and SetUpFirstExpList().

◆ m_inputfiles

std::map<std::string, std::vector<std::string> > Nektar::FieldUtils::Field::m_inputfiles

Definition at line 102 of file Field.hpp.

◆ m_nParts

int Nektar::FieldUtils::Field::m_nParts = 1

Definition at line 97 of file Field.hpp.

◆ m_numHomogeneousDir

int Nektar::FieldUtils::Field::m_numHomogeneousDir

Definition at line 85 of file Field.hpp.

Referenced by CreateExp().

◆ m_partComm

LibUtilities::CommSharedPtr Nektar::FieldUtils::Field::m_partComm

Definition at line 96 of file Field.hpp.

◆ m_requireBoundaryExpansion

bool Nektar::FieldUtils::Field::m_requireBoundaryExpansion

Definition at line 90 of file Field.hpp.

Referenced by AppendExpList().

◆ m_session

LibUtilities::SessionReaderSharedPtr Nektar::FieldUtils::Field::m_session

Definition at line 100 of file Field.hpp.

Referenced by AppendExpList(), ClearField(), CreateExp(), and SetUpFirstExpList().

◆ m_useFFT

bool Nektar::FieldUtils::Field::m_useFFT

Definition at line 92 of file Field.hpp.

Referenced by AppendExpList(), and SetUpFirstExpList().

◆ m_variables

std::vector<std::string> Nektar::FieldUtils::Field::m_variables

Definition at line 83 of file Field.hpp.

Referenced by ClearField(), and CreateExp().

◆ m_verbose

bool Nektar::FieldUtils::Field::m_verbose

Definition at line 79 of file Field.hpp.

Referenced by CreateExp().

◆ m_vm

po::variables_map Nektar::FieldUtils::Field::m_vm

Definition at line 98 of file Field.hpp.

◆ m_writeBndFld

bool Nektar::FieldUtils::Field::m_writeBndFld

Definition at line 104 of file Field.hpp.