Nektar++
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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:183
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 680 of file Field.hpp.

684  {
685  if (var.compare("DefaultVar") == 0 && m_requireBoundaryExpansion)
686  {
687  if (m_session->GetVariables().size())
688  {
689  var = m_session->GetVariables()[0];
690  }
691  }
693  switch (m_graph->GetMeshDimension())
694  {
695  case 1:
696  {
697  if (NumHomogeneousDir == 1)
698  {
701  "ContFieldHomogeneous1D or "
702  "DisContFieldHomogenenous1D has not been "
703  "implemented");
704 
706  std::dynamic_pointer_cast<
707  MultiRegions::ExpList2DHomogeneous1D>(m_exp[0]);
708 
710  AllocateSharedPtr(*tmp2);
711  }
712  else if (NumHomogeneousDir == 2)
713  {
715  {
717  std::dynamic_pointer_cast<
718  MultiRegions::ContField3DHomogeneous2D>(
719  m_exp[0]);
720 
721  tmp = MemoryManager<
722  MultiRegions::ContField3DHomogeneous2D>::
723  AllocateSharedPtr(*tmp2);
724  }
726  {
728  tmp2 = std::dynamic_pointer_cast<
729  MultiRegions::DisContField3DHomogeneous2D>(
730  m_exp[0]);
731 
732  tmp = MemoryManager<
733  MultiRegions::DisContField3DHomogeneous2D>::
734  AllocateSharedPtr(*tmp2);
735  }
736  else
737  {
739  std::dynamic_pointer_cast<
740  MultiRegions::ExpList3DHomogeneous2D>(m_exp[0]);
741 
742  tmp = MemoryManager<
743  MultiRegions::ExpList3DHomogeneous2D>::
744  AllocateSharedPtr(*tmp2);
745  }
746  }
747  else
748  {
750  {
752  std::dynamic_pointer_cast<
753  MultiRegions::ContField>(m_exp[0]);
754 
757  }
759  {
761  std::dynamic_pointer_cast<
762  MultiRegions::DisContField>(m_exp[0]);
763 
766  }
767  else
768  {
770  std::dynamic_pointer_cast<
771  MultiRegions::ExpList>(m_exp[0]);
772 
773  tmp = MemoryManager<
774  MultiRegions::ExpList>::AllocateSharedPtr(*tmp2);
775  }
776  }
777  }
778  break;
779  case 2:
780  {
781  if (NumHomogeneousDir == 1)
782  {
784  {
785  if (NewField)
786  {
787  bool dealiasing = false;
788 
789  tmp = MemoryManager<
790  MultiRegions::ContField3DHomogeneous1D>::
791  AllocateSharedPtr(
792  m_session, m_exp[0]
793  ->GetHomogeneousBasis()
794  ->GetBasisKey(),
795  m_exp[0]->GetHomoLen(), m_useFFT,
796  dealiasing, m_graph, var);
797  }
798  else
799  {
801  tmp2 = std::dynamic_pointer_cast<
802  MultiRegions::ContField3DHomogeneous1D>(
803  m_exp[0]);
804 
805  ASSERTL0(tmp2, "Failed to type cast m_exp[0]");
806  tmp = MemoryManager<
807  MultiRegions::ContField3DHomogeneous1D>::
808  AllocateSharedPtr(*tmp2, m_graph, var);
809  }
810  }
812  {
813  if (NewField)
814  {
815  bool dealiasing = false;
816 
817  tmp = MemoryManager<
818  MultiRegions::DisContField3DHomogeneous1D>::
819  AllocateSharedPtr(
820  m_session, m_exp[0]
821  ->GetHomogeneousBasis()
822  ->GetBasisKey(),
823  m_exp[0]->GetHomoLen(), m_useFFT,
824  dealiasing, m_graph, var);
825  }
826  else
827  {
829  tmp2 = std::dynamic_pointer_cast<
830  MultiRegions::DisContField3DHomogeneous1D>(
831  m_exp[0]);
832  ASSERTL0(tmp2, "Failed to type cast m_exp[0]");
833 
834  tmp = MemoryManager<
835  MultiRegions::DisContField3DHomogeneous1D>::
836  AllocateSharedPtr(*tmp2);
837  }
838  }
839  else
840  {
841  if (NewField)
842  {
843  bool dealiasing = false;
844 
845  tmp = MemoryManager<
846  MultiRegions::ExpList3DHomogeneous1D>::
847  AllocateSharedPtr(
848  m_session, m_exp[0]
849  ->GetHomogeneousBasis()
850  ->GetBasisKey(),
851  m_exp[0]->GetHomoLen(), m_useFFT,
852  dealiasing, m_graph);
853  }
854  else
855  {
857  std::dynamic_pointer_cast<
858  MultiRegions::ExpList3DHomogeneous1D>(
859  m_exp[0]);
860  ASSERTL0(tmp2, "Failed to type cast m_exp[0]");
861 
862  tmp = MemoryManager<
863  MultiRegions::ExpList3DHomogeneous1D>::
864  AllocateSharedPtr(*tmp2);
865  }
866  }
867  }
868  else
869  {
871  {
872  if (NewField)
873  {
876  }
877  else // call copy constructor
878  {
879 
881  std::dynamic_pointer_cast<
882  MultiRegions::ContField>(m_exp[0]);
883 
885  AllocateSharedPtr(*tmp2, m_graph, var);
886  }
887  }
889  {
890  if (NewField)
891  {
894  }
895  else // call copy constructor
896  {
898  std::dynamic_pointer_cast<
899  MultiRegions::DisContField>(m_exp[0]);
900 
902  AllocateSharedPtr(*tmp2, m_graph, var);
903  }
904  }
905  else
906  {
908  std::dynamic_pointer_cast<
909  MultiRegions::ExpList>(m_exp[0]);
910 
911  tmp = MemoryManager<
912  MultiRegions::ExpList>::AllocateSharedPtr(*tmp2);
913  }
914  }
915  }
916  break;
917  case 3:
918  {
920  {
921  if (NewField)
922  {
925  }
926  else
927  {
929  std::dynamic_pointer_cast<
930  MultiRegions::ContField>(m_exp[0]);
931 
933  AllocateSharedPtr(*tmp2, m_graph, var);
934  }
935  }
937  {
938  if (NewField)
939  {
942  }
943  else
944  {
946  std::dynamic_pointer_cast<
947  MultiRegions::DisContField>(m_exp[0]);
948 
950  AllocateSharedPtr(*tmp2, m_graph, var);
951  }
952  }
953  else
954  {
956  std::dynamic_pointer_cast<MultiRegions::ExpList>(
957  m_exp[0]);
958 
959  tmp = MemoryManager<
960  MultiRegions::ExpList>::AllocateSharedPtr(*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:216
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:385
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:292
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, 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:174
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.

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

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 659 of file Field.hpp.

661  {
663  std::string fmt = LibUtilities::FieldIO::GetFileType(filename, c);
664 
665  auto it = m_fld.find(fmt);
666 
667  if (it == m_fld.end())
668  {
671  m_fld[fmt] = fld;
672  return fld;
673  }
674  else
675  {
676  return it->second;
677  }
678  }
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:145
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:306
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 328 of file Field.hpp.

330  {
331 
333 
334  // Set up expansion list
335  int expdim = m_graph->GetMeshDimension();
336  bool dealiasing = false;
337 
338  m_session->MatchSolverInfo("USEFFT", "FFTW", m_useFFT, false);
339 
340  switch (expdim)
341  {
342  case 1:
343  {
344  ASSERTL0(NumHomogeneousDir <= 2,
345  "Quasi-3D approach is only set up for 1 or 2 "
346  "homogeneous directions");
347 
348  if (NumHomogeneousDir == 1)
349  {
351 
352  // Define Homogeneous expansion
353  int nplanes;
354  NekDouble ly;
356 
357  if (fldfilegiven)
358  {
359  nplanes = m_fielddef[0]->m_numModes[1];
360  ly = m_fielddef[0]->m_homogeneousLengths[0];
361  btype = m_fielddef[0]->m_basis[1];
362  }
363  else
364  {
365  m_session->LoadParameter("HomModesZ", nplanes);
366  m_session->LoadParameter("LY", ly);
367  btype = LibUtilities::eFourier;
368  }
369 
370  // Choose points to be at evenly spaced points at
371  // nplanes points
372  const LibUtilities::PointsKey Pkey(
374 
375  const LibUtilities::BasisKey Bkey(btype, nplanes, Pkey);
376 
379  {
380  ASSERTL0(false, "ContFieldHomogeneous1D or "
381  "DisContFieldHomogenenous1D has "
382  "not been implemented");
383  }
384 
385  Exp2DH1 =
388  dealiasing, m_graph,
390  exp = Exp2DH1;
391  }
392  else if (NumHomogeneousDir == 2)
393  {
395 
396  int nylines, nzlines;
397  NekDouble ly, lz;
398  LibUtilities::BasisType btype1, btype2;
399 
400  if (fldfilegiven)
401  {
402  nylines = m_fielddef[0]->m_numModes[1];
403  nzlines = m_fielddef[0]->m_numModes[2];
404  ly = m_fielddef[0]->m_homogeneousLengths[0];
405  lz = m_fielddef[0]->m_homogeneousLengths[1];
406  btype1 = m_fielddef[0]->m_basis[1];
407  btype2 = m_fielddef[0]->m_basis[2];
408  }
409  else
410  {
411  m_session->LoadParameter("HomModesY", nylines);
412  m_session->LoadParameter("HomModesZ", nzlines);
413  m_session->LoadParameter("LY", ly);
414  m_session->LoadParameter("LZ", lz);
415  btype1 = LibUtilities::eFourier;
416  btype2 = LibUtilities::eFourier;
417  }
418 
419  // Choose points to be at evenly spaced points at
420  // nplanes points
421  const LibUtilities::PointsKey PkeyY(
423  const LibUtilities::BasisKey BkeyY(btype1, nylines, PkeyY);
424 
425  const LibUtilities::PointsKey PkeyZ(
427  const LibUtilities::BasisKey BkeyZ(btype2, nzlines, PkeyZ);
428 
430  {
431  Exp3DH2 = MemoryManager<
432  MultiRegions::ContField3DHomogeneous2D>::
433  AllocateSharedPtr(m_session, BkeyY, BkeyZ, ly, lz,
434  m_useFFT, dealiasing, m_graph,
435  m_session->GetVariable(0),
437  }
439  {
440  Exp3DH2 = MemoryManager<
441  MultiRegions::DisContField3DHomogeneous2D>::
442  AllocateSharedPtr(m_session, BkeyY, BkeyZ, ly, lz,
443  m_useFFT, dealiasing, m_graph,
444  m_session->GetVariable(0),
446  }
447  else
448  {
449  Exp3DH2 = MemoryManager<
450  MultiRegions::ExpList3DHomogeneous2D>::
451  AllocateSharedPtr(m_session, BkeyY, BkeyZ, ly, lz,
452  m_useFFT, dealiasing, m_graph,
454  }
455 
456  exp = Exp3DH2;
457  }
458  else
459  {
461 
463  {
466  m_session->GetVariable(0),
467  true,false,
469  }
471  {
474  m_session->GetVariable(0),
475  true,
477  }
478  else
479  {
482  true,"DefaultVar",
484  }
485 
486  exp = Exp1D;
487  }
488  }
489  break;
490  case 2:
491  {
492  ASSERTL0(NumHomogeneousDir <= 1,
493  "NumHomogeneousDir is only set up for 1");
494 
495  if (NumHomogeneousDir == 1)
496  {
498 
499  // Define Homogeneous expansion
500  int nplanes;
501  NekDouble lz;
505 
506  if (fldfilegiven)
507  {
508  nplanes = m_fielddef[0]->m_numModes[2];
509  lz = m_fielddef[0]->m_homogeneousLengths[0];
510  btype = m_fielddef[0]->m_basis[2];
511 
513  {
514  btype = LibUtilities::eFourier;
515  m_fielddef[0]->m_basis[2] =
517  if (nplanes <= 2)
518  {
519  nplanes = 4;
520  }
521  }
522  else if (btype == LibUtilities::eFourierHalfModeRe &&
523  nplanes == 1)
524  {
526  }
527  }
528  else
529  {
530  m_session->LoadParameter("HomModesZ", nplanes);
531  m_session->LoadParameter("LZ", lz);
532  btype = LibUtilities::eFourier;
533  }
534  // Choose points to be at evenly spaced points at
535  // nplanes points
536  const LibUtilities::PointsKey Pkey(nplanes, ptype);
537 
538  const LibUtilities::BasisKey Bkey(btype, nplanes, Pkey);
539 
541  {
542  Exp3DH1 = MemoryManager<
543  MultiRegions::ContField3DHomogeneous1D>::
544  AllocateSharedPtr(m_session, Bkey, lz, m_useFFT,
545  dealiasing, m_graph,
546  m_session->GetVariable(0),
548  }
549 
551  {
552  Exp3DH1 = MemoryManager<
553  MultiRegions::DisContField3DHomogeneous1D>::
554  AllocateSharedPtr(m_session, Bkey, lz, m_useFFT,
555  dealiasing, m_graph,
556  m_session->GetVariable(0),
558  }
559 
560  else
561  {
562  Exp3DH1 = MemoryManager<
563  MultiRegions::ExpList3DHomogeneous1D>::
564  AllocateSharedPtr(m_session, Bkey, lz, m_useFFT,
565  dealiasing, m_graph,
566  "DefaultVar",
568  }
569  exp = Exp3DH1;
570  }
571  else
572  {
574 
576  {
579  m_session->GetVariable(0),
580  true,false,
582  }
583 
585  {
588  m_session->GetVariable(0),
589  true,true,
591  }
592  else
593  {
596  true,
597  "DefaultVar",
599  }
600 
601  exp = Exp2D;
602  }
603  }
604  break;
605  case 3:
606  {
608 
610  {
613  m_session->GetVariable(0),
614  true,
615  false,
617  }
619  {
622  m_session->GetVariable(0),
623  true,
625  }
626  else
627  {
628  Exp3D = MemoryManager<
629  MultiRegions::ExpList>::AllocateSharedPtr(
630  m_session,
631  m_graph,
632  true,
633  "DefaultVar",
635  }
636 
637  exp = Exp3D;
638  }
639  break;
640  default:
641  ASSERTL0(false, "Expansion dimension not recognised");
642  break;
643  }
644 
645  return exp;
646  };
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:65
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:64
@ eFourierSingleMode
Fourier ModifiedExpansion with just the first mode .
Definition: BasisType.h:59
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode
Definition: BasisType.h:60
@ eFourier
Fourier Expansion .
Definition: BasisType.h:53

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.