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

Constructor & Destructor Documentation

◆ Field()

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

Definition at line 62 of file Field.hpp.

67 {
68 }
static PtsFieldSharedPtr NullPtsField
Definition: PtsField.h:185
bool m_declareExpansionAsDisContField
Definition: Field.hpp:87
bool m_declareExpansionAsContField
Definition: Field.hpp:86
LibUtilities::PtsFieldSharedPtr m_fieldPts
Definition: Field.hpp:107

◆ ~Field()

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

Definition at line 70 of file Field.hpp.

71 {
72 if (m_comm)
73 {
74 m_comm->Finalise();
75 }
76 }
LibUtilities::CommSharedPtr m_comm
Definition: Field.hpp:93

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

666 {
667 if (var.compare("DefaultVar") == 0 && m_requireBoundaryExpansion)
668 {
669 if (m_session->GetVariables().size())
670 {
671 var = m_session->GetVariables()[0];
672 }
673 }
675 switch (m_graph->GetMeshDimension())
676 {
677 case 1:
678 {
679 if (NumHomogeneousDir == 1)
680 {
683 "ContFieldHomogeneous1D or "
684 "DisContFieldHomogenenous1D has not been "
685 "implemented");
686
688 std::dynamic_pointer_cast<
689 MultiRegions::ExpList2DHomogeneous1D>(m_exp[0]);
690
692 AllocateSharedPtr(*tmp2);
693 }
694 else if (NumHomogeneousDir == 2)
695 {
697 {
699 std::dynamic_pointer_cast<
700 MultiRegions::ContField3DHomogeneous2D>(
701 m_exp[0]);
702
703 tmp = MemoryManager<
704 MultiRegions::ContField3DHomogeneous2D>::
705 AllocateSharedPtr(*tmp2);
706 }
708 {
710 tmp2 = std::dynamic_pointer_cast<
711 MultiRegions::DisContField3DHomogeneous2D>(
712 m_exp[0]);
713
714 tmp = MemoryManager<
715 MultiRegions::DisContField3DHomogeneous2D>::
716 AllocateSharedPtr(*tmp2);
717 }
718 else
719 {
721 std::dynamic_pointer_cast<
722 MultiRegions::ExpList3DHomogeneous2D>(m_exp[0]);
723
724 tmp = MemoryManager<
725 MultiRegions::ExpList3DHomogeneous2D>::
726 AllocateSharedPtr(*tmp2);
727 }
728 }
729 else
730 {
732 {
734 std::dynamic_pointer_cast<MultiRegions::ContField>(
735 m_exp[0]);
736
739 false,
741 }
743 {
745 std::dynamic_pointer_cast<
746 MultiRegions::DisContField>(m_exp[0]);
747
751 }
752 else
753 {
755 std::dynamic_pointer_cast<MultiRegions::ExpList>(
756 m_exp[0]);
757
758 tmp = MemoryManager<
759 MultiRegions::ExpList>::AllocateSharedPtr(*tmp2);
760 }
761 }
762 }
763 break;
764 case 2:
765 {
766 if (NumHomogeneousDir == 1)
767 {
769 {
770 if (NewField)
771 {
772 bool dealiasing = false;
773
774 tmp = MemoryManager<
775 MultiRegions::ContField3DHomogeneous1D>::
776 AllocateSharedPtr(m_session,
777 m_exp[0]
778 ->GetHomogeneousBasis()
779 ->GetBasisKey(),
780 m_exp[0]->GetHomoLen(),
781 m_useFFT, dealiasing, m_graph,
782 var,
784 }
785 else
786 {
788 tmp2 = std::dynamic_pointer_cast<
789 MultiRegions::ContField3DHomogeneous1D>(
790 m_exp[0]);
791
792 ASSERTL0(tmp2, "Failed to type cast m_exp[0]");
793 tmp = MemoryManager<
794 MultiRegions::ContField3DHomogeneous1D>::
795 AllocateSharedPtr(*tmp2, m_graph, var);
796 }
797 }
799 {
800 if (NewField)
801 {
802 bool dealiasing = false;
803
804 tmp = MemoryManager<
805 MultiRegions::DisContField3DHomogeneous1D>::
806 AllocateSharedPtr(m_session,
807 m_exp[0]
808 ->GetHomogeneousBasis()
809 ->GetBasisKey(),
810 m_exp[0]->GetHomoLen(),
811 m_useFFT, dealiasing, m_graph,
812 var,
814 }
815 else
816 {
818 tmp2 = std::dynamic_pointer_cast<
819 MultiRegions::DisContField3DHomogeneous1D>(
820 m_exp[0]);
821 ASSERTL0(tmp2, "Failed to type cast m_exp[0]");
822
823 tmp = MemoryManager<
824 MultiRegions::DisContField3DHomogeneous1D>::
825 AllocateSharedPtr(*tmp2);
826 }
827 }
828 else
829 {
830 if (NewField)
831 {
832 bool dealiasing = false;
833
834 tmp = MemoryManager<
835 MultiRegions::ExpList3DHomogeneous1D>::
836 AllocateSharedPtr(m_session,
837 m_exp[0]
838 ->GetHomogeneousBasis()
839 ->GetBasisKey(),
840 m_exp[0]->GetHomoLen(),
841 m_useFFT, dealiasing, m_graph,
842 "DefaultVar",
844 }
845 else
846 {
848 std::dynamic_pointer_cast<
849 MultiRegions::ExpList3DHomogeneous1D>(
850 m_exp[0]);
851 ASSERTL0(tmp2, "Failed to type cast m_exp[0]");
852
853 tmp = MemoryManager<
854 MultiRegions::ExpList3DHomogeneous1D>::
855 AllocateSharedPtr(*tmp2);
856 }
857 }
858 }
859 else
860 {
862 {
863 if (NewField)
864 {
867 false,
869 }
870 else // call copy constructor
871 {
872
874 std::dynamic_pointer_cast<
875 MultiRegions::ContField>(m_exp[0]);
876
878 AllocateSharedPtr(*tmp2, m_graph, var);
879 }
880 }
882 {
883 if (NewField)
884 {
887 true,
889 }
890 else // call copy constructor
891 {
893 std::dynamic_pointer_cast<
894 MultiRegions::DisContField>(m_exp[0]);
895
897 AllocateSharedPtr(*tmp2, m_graph, var);
898 }
899 }
900 else
901 {
903 std::dynamic_pointer_cast<MultiRegions::ExpList>(
904 m_exp[0]);
905
906 tmp = MemoryManager<
907 MultiRegions::ExpList>::AllocateSharedPtr(*tmp2);
908 }
909 }
910 }
911 break;
912 case 3:
913 {
915 {
916 if (NewField)
917 {
920 false,
922 }
923 else
924 {
926 std::dynamic_pointer_cast<MultiRegions::ContField>(
927 m_exp[0]);
928
929 tmp = MemoryManager<
930 MultiRegions::ContField>::AllocateSharedPtr(*tmp2,
931 m_graph,
932 var);
933 }
934 }
936 {
937 if (NewField)
938 {
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 =
961 *tmp2);
962 }
963 }
964 break;
965 default:
966 ASSERTL0(false, "Expansion dimension not recognised");
967 break;
968 }
969
970 return tmp;
971 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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:345
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:268
SpatialDomains::MeshGraphSharedPtr m_graph
Definition: Field.hpp:100
LibUtilities::SessionReaderSharedPtr m_session
Definition: Field.hpp:99
std::vector< MultiRegions::ExpListSharedPtr > m_exp
Definition: Field.hpp:80

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

974 {
978 m_exp.clear();
979 m_fielddef = std::vector<LibUtilities::FieldDefinitionsSharedPtr>();
980 m_data = std::vector<std::vector<NekDouble>>();
981 m_variables.clear();
982 }
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::vector< std::string > m_variables
Definition: Field.hpp:82
std::vector< LibUtilities::FieldDefinitionsSharedPtr > m_fielddef
Definition: Field.hpp:78
std::vector< std::vector< double > > m_data
Definition: Field.hpp:79

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

◆ CreateExp()

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

Definition at line 123 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 =
142 fldfilegiven && !vm.count("use-session-expansion");
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 (size_t 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
198 return;
199 }
200
201 // Adjust number of quadrature points
202 if (vm.count("output-points"))
203 {
204 int nPointsNew = vm["output-points"].as<int>();
205 m_graph->SetExpansionInfoToPointOrder(nPointsNew);
206 }
207
208 if (m_verbose)
209 {
210 if (m_comm->TreatAsRankZero())
211 {
212 timerpart.Stop();
213 NekDouble cpuTime = timerpart.TimePerTest(1);
214
215 std::stringstream ss;
216 ss << cpuTime << "s";
217 std::cout << "\t CreateExp setexpansion CPU Time: "
218 << std::setw(8) << std::left << ss.str()
219 << std::endl;
220 timerpart.Start();
221 }
222 }
223 // Override number of planes with value from cmd line
224 if (m_numHomogeneousDir == 1 && vm.count("output-points-hom-z"))
225 {
226 int expdim = m_graph->GetMeshDimension();
227 m_fielddef[0]->m_numModes[expdim] =
228 vm["output-points-hom-z"].as<int>();
229 }
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 << ss1.str()
243 << std::endl;
244 }
245 }
246 }
247
248 if (fldfilegiven)
249 {
250 size_t 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("use-session-variables"))
255 {
256 m_variables = vars;
257 }
258 nfields = m_variables.size();
259
260 m_exp.resize(nfields * nstrips);
261 // declare other fields;
262 for (size_t 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] =
273 }
274 }
275 else
276 {
277 if (vars.size())
278 {
279 m_exp[s * nfields + i] =
281 }
282 else
283 {
284 m_exp[s * nfields + i] =
286 }
287 }
288 }
289 }
290
291 // Extract data to coeffs and bwd transform
292 for (size_t s = 0; s < nstrips; ++s) // homogeneous strip varient
293 {
294 for (j = 0; j < nfields; ++j)
295 {
296 for (i = 0; i < m_data.size() / nstrips; ++i)
297 {
298 size_t n = i * nstrips + s;
299 // In case of multiple flds, we might not have a
300 // variable in this m_data[n] -> skip in this case
301 auto it =
302 find(m_fielddef[n]->m_fields.begin(),
303 m_fielddef[n]->m_fields.end(), m_variables[j]);
304 if (it != m_fielddef[n]->m_fields.end())
305 {
306 m_exp[s * nfields + j]->ExtractDataToCoeffs(
307 m_fielddef[n], m_data[n], m_variables[j],
308 m_exp[s * nfields + j]->UpdateCoeffs());
309 }
310 }
311 m_exp[s * nfields + j]->BwdTrans(
312 m_exp[s * nfields + j]->GetCoeffs(),
313 m_exp[s * nfields + j]->UpdatePhys());
314 }
315 }
316 // Clear fielddef and data
317 // (they should not be used after running this module)
318 m_fielddef = std::vector<LibUtilities::FieldDefinitionsSharedPtr>();
319 m_data = std::vector<std::vector<NekDouble>>();
320 }
321 }
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:447
std::vector< double > d(NPUPPER *NPUPPER)
double NekDouble
FIELD_UTILS_EXPORT MultiRegions::ExpListSharedPtr SetUpFirstExpList(int NumHomogeneousDir, bool fldfilegiven=false)
Definition: Field.hpp:323
FIELD_UTILS_EXPORT MultiRegions::ExpListSharedPtr AppendExpList(int NumHomogeneousDir, std::string var="DefaultVar", bool NewField=false)
Definition: Field.hpp:663

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), AppendExpList(), Nektar::UnitTests::d(), 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 642 of file Field.hpp.

644 {
646 std::string fmt = LibUtilities::FieldIO::GetFileType(filename, c);
647
648 auto it = m_fld.find(fmt);
649
650 if (it == m_fld.end())
651 {
654 m_fld[fmt] = fld;
655 return fld;
656 }
657 else
658 {
659 return it->second;
660 }
661 }
static const std::string GetFileType(const std::string &filename, CommSharedPtr comm)
Determine file type of given input file.
Definition: FieldIO.cpp:95
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:322
FieldIOFactory & GetFieldIOFactory()
Returns the FieldIO factory.
Definition: FieldIO.cpp:70
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
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:987

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

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

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

325 {
326
328
329 // Set up expansion list
330 int expdim = m_graph->GetMeshDimension();
331 bool dealiasing = false;
332
333 m_session->MatchSolverInfo("USEFFT", "FFTW", m_useFFT, false);
334
335 switch (expdim)
336 {
337 case 1:
338 {
339 ASSERTL0(NumHomogeneousDir <= 2,
340 "Quasi-3D approach is only set up for 1 or 2 "
341 "homogeneous directions");
342
343 if (NumHomogeneousDir == 1)
344 {
346
347 // Define Homogeneous expansion
348 int nplanes;
349 NekDouble ly;
351
352 if (fldfilegiven)
353 {
354 nplanes = m_fielddef[0]->m_numModes[1];
355 ly = m_fielddef[0]->m_homogeneousLengths[0];
356 btype = m_fielddef[0]->m_basis[1];
357 }
358 else
359 {
360 m_session->LoadParameter("HomModesZ", nplanes);
361 m_session->LoadParameter("LY", ly);
363 }
364
365 // Choose points to be at evenly spaced points at
366 // nplanes points
367 const LibUtilities::PointsKey Pkey(
369
370 const LibUtilities::BasisKey Bkey(btype, nplanes, Pkey);
371
374 {
375 ASSERTL0(false, "ContFieldHomogeneous1D or "
376 "DisContFieldHomogenenous1D has "
377 "not been implemented");
378 }
379
380 Exp2DH1 =
383 dealiasing, m_graph,
385 exp = Exp2DH1;
386 }
387 else if (NumHomogeneousDir == 2)
388 {
390
391 int nylines, nzlines;
392 NekDouble ly, lz;
393 LibUtilities::BasisType btype1, btype2;
394
395 if (fldfilegiven)
396 {
397 nylines = m_fielddef[0]->m_numModes[1];
398 nzlines = m_fielddef[0]->m_numModes[2];
399 ly = m_fielddef[0]->m_homogeneousLengths[0];
400 lz = m_fielddef[0]->m_homogeneousLengths[1];
401 btype1 = m_fielddef[0]->m_basis[1];
402 btype2 = m_fielddef[0]->m_basis[2];
403 }
404 else
405 {
406 m_session->LoadParameter("HomModesY", nylines);
407 m_session->LoadParameter("HomModesZ", nzlines);
408 m_session->LoadParameter("LY", ly);
409 m_session->LoadParameter("LZ", lz);
410 btype1 = LibUtilities::eFourier;
411 btype2 = LibUtilities::eFourier;
412 }
413
414 // Choose points to be at evenly spaced points at
415 // nplanes points
416 const LibUtilities::PointsKey PkeyY(
418 const LibUtilities::BasisKey BkeyY(btype1, nylines, PkeyY);
419
420 const LibUtilities::PointsKey PkeyZ(
422 const LibUtilities::BasisKey BkeyZ(btype2, nzlines, PkeyZ);
423
425 {
426 Exp3DH2 = MemoryManager<
427 MultiRegions::ContField3DHomogeneous2D>::
428 AllocateSharedPtr(m_session, BkeyY, BkeyZ, ly, lz,
429 m_useFFT, dealiasing, m_graph,
430 m_session->GetVariable(0),
432 }
434 {
435 Exp3DH2 = MemoryManager<
436 MultiRegions::DisContField3DHomogeneous2D>::
437 AllocateSharedPtr(m_session, BkeyY, BkeyZ, ly, lz,
438 m_useFFT, dealiasing, m_graph,
439 m_session->GetVariable(0),
441 }
442 else
443 {
444 Exp3DH2 = MemoryManager<
445 MultiRegions::ExpList3DHomogeneous2D>::
446 AllocateSharedPtr(m_session, BkeyY, BkeyZ, ly, lz,
447 m_useFFT, dealiasing, m_graph,
449 }
450
451 exp = Exp3DH2;
452 }
453 else
454 {
456
458 {
461 m_session, m_graph, m_session->GetVariable(0),
462 true, false, Collections::eNoCollection);
463 }
465 {
468 m_session->GetVariable(0), true,
470 }
471 else
472 {
475 "DefaultVar",
477 }
478
479 exp = Exp1D;
480 }
481 }
482 break;
483 case 2:
484 {
485 ASSERTL0(NumHomogeneousDir <= 1,
486 "NumHomogeneousDir is only set up for 1");
487
488 if (NumHomogeneousDir == 1)
489 {
491
492 // Define Homogeneous expansion
493 int nplanes;
494 NekDouble lz;
498
499 if (fldfilegiven)
500 {
501 nplanes = m_fielddef[0]->m_numModes[2];
502 lz = m_fielddef[0]->m_homogeneousLengths[0];
503 btype = m_fielddef[0]->m_basis[2];
504
506 {
508 m_fielddef[0]->m_basis[2] =
510 if (nplanes <= 2)
511 {
512 nplanes = 4;
513 }
514 }
515 else if (btype == LibUtilities::eFourierHalfModeRe &&
516 nplanes == 1)
517 {
519 }
520 }
521 else
522 {
523 m_session->LoadParameter("HomModesZ", nplanes);
524 m_session->LoadParameter("LZ", lz);
526 }
527 // Choose points to be at evenly spaced points at
528 // nplanes points
529 const LibUtilities::PointsKey Pkey(nplanes, ptype);
530
531 const LibUtilities::BasisKey Bkey(btype, nplanes, Pkey);
532
534 {
535 Exp3DH1 = MemoryManager<
536 MultiRegions::ContField3DHomogeneous1D>::
537 AllocateSharedPtr(m_session, Bkey, lz, m_useFFT,
538 dealiasing, m_graph,
539 m_session->GetVariable(0),
541 }
542
544 {
545 Exp3DH1 = MemoryManager<
546 MultiRegions::DisContField3DHomogeneous1D>::
547 AllocateSharedPtr(m_session, Bkey, lz, m_useFFT,
548 dealiasing, m_graph,
549 m_session->GetVariable(0),
551 }
552
553 else
554 {
555 Exp3DH1 = MemoryManager<
556 MultiRegions::ExpList3DHomogeneous1D>::
557 AllocateSharedPtr(m_session, Bkey, lz, m_useFFT,
558 dealiasing, m_graph, "DefaultVar",
560 }
561 exp = Exp3DH1;
562 }
563 else
564 {
566
568 {
571 m_session, m_graph, m_session->GetVariable(0),
572 true, false, Collections::eNoCollection);
573 }
574
576 {
579 m_session->GetVariable(0), true,
581 }
582 else
583 {
586 "DefaultVar",
588 }
589
590 exp = Exp2D;
591 }
592 }
593 break;
594 case 3:
595 {
597
599 {
602 m_session->GetVariable(0), true,
604 }
606 {
609 m_session->GetVariable(0), true,
611 }
612 else
613 {
614 Exp3D =
616 m_session, m_graph, true, "DefaultVar",
618 }
619
620 exp = Exp3D;
621 }
622 break;
623 default:
624 ASSERTL0(false, "Expansion dimension not recognised");
625 break;
626 }
627
628 return exp;
629 };
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:74
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:73
@ eFourierSingleMode
Fourier ModifiedExpansion with just the first mode .
Definition: BasisType.h:64
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode .
Definition: BasisType.h:66
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55

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

◆ m_bndRegionsToWrite

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

Definition at line 104 of file Field.hpp.

◆ m_comm

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

Definition at line 93 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 79 of file Field.hpp.

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

◆ m_declareExpansionAsContField

bool Nektar::FieldUtils::Field::m_declareExpansionAsContField

Definition at line 86 of file Field.hpp.

Referenced by AppendExpList(), and SetUpFirstExpList().

◆ m_declareExpansionAsDisContField

bool Nektar::FieldUtils::Field::m_declareExpansionAsDisContField

Definition at line 87 of file Field.hpp.

Referenced by AppendExpList(), and SetUpFirstExpList().

◆ m_defComm

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

Definition at line 94 of file Field.hpp.

◆ m_exp

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

Definition at line 80 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 78 of file Field.hpp.

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

◆ m_fieldMetaDataMap

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

Definition at line 109 of file Field.hpp.

◆ m_fieldPts

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

Definition at line 107 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 987 of file Field.hpp.

Referenced by FieldIOForFile().

◆ m_graph

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

Definition at line 100 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 101 of file Field.hpp.

◆ m_nParts

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

Definition at line 96 of file Field.hpp.

◆ m_numHomogeneousDir

int Nektar::FieldUtils::Field::m_numHomogeneousDir

Definition at line 84 of file Field.hpp.

Referenced by CreateExp().

◆ m_partComm

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

Definition at line 95 of file Field.hpp.

◆ m_requireBoundaryExpansion

bool Nektar::FieldUtils::Field::m_requireBoundaryExpansion

Definition at line 89 of file Field.hpp.

Referenced by AppendExpList().

◆ m_session

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

Definition at line 99 of file Field.hpp.

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

◆ m_useFFT

bool Nektar::FieldUtils::Field::m_useFFT

Definition at line 91 of file Field.hpp.

Referenced by AppendExpList(), and SetUpFirstExpList().

◆ m_variables

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

Definition at line 82 of file Field.hpp.

Referenced by ClearField(), and CreateExp().

◆ m_verbose

bool Nektar::FieldUtils::Field::m_verbose

Definition at line 77 of file Field.hpp.

Referenced by CreateExp().

◆ m_vm

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

Definition at line 97 of file Field.hpp.

◆ m_writeBndFld

bool Nektar::FieldUtils::Field::m_writeBndFld

Definition at line 103 of file Field.hpp.