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:89
bool m_declareExpansionAsContField
Definition: Field.hpp:88
LibUtilities::PtsFieldSharedPtr m_fieldPts
Definition: Field.hpp:109

◆ ~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:95

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

667 {
668 if (var.compare("DefaultVar") == 0 && m_requireBoundaryExpansion)
669 {
670 if (m_session->GetVariables().size())
671 {
672 var = m_session->GetVariables()[0];
673 }
674 }
676 switch (m_graph->GetMeshDimension())
677 {
678 case 1:
679 {
680 if (NumHomogeneousDir == 1)
681 {
684 "ContFieldHomogeneous1D or "
685 "DisContFieldHomogenenous1D has not been "
686 "implemented");
687
689 std::dynamic_pointer_cast<
690 MultiRegions::ExpList2DHomogeneous1D>(m_exp[0]);
691
693 AllocateSharedPtr(*tmp2);
694 }
695 else if (NumHomogeneousDir == 2)
696 {
698 {
700 std::dynamic_pointer_cast<
701 MultiRegions::ContField3DHomogeneous2D>(
702 m_exp[0]);
703
704 tmp = MemoryManager<
705 MultiRegions::ContField3DHomogeneous2D>::
706 AllocateSharedPtr(*tmp2);
707 }
709 {
711 tmp2 = std::dynamic_pointer_cast<
712 MultiRegions::DisContField3DHomogeneous2D>(
713 m_exp[0]);
714
715 tmp = MemoryManager<
716 MultiRegions::DisContField3DHomogeneous2D>::
717 AllocateSharedPtr(*tmp2);
718 }
719 else
720 {
722 std::dynamic_pointer_cast<
723 MultiRegions::ExpList3DHomogeneous2D>(m_exp[0]);
724
725 tmp = MemoryManager<
726 MultiRegions::ExpList3DHomogeneous2D>::
727 AllocateSharedPtr(*tmp2);
728 }
729 }
730 else
731 {
733 {
735 std::dynamic_pointer_cast<MultiRegions::ContField>(
736 m_exp[0]);
737
740 false,
742 }
744 {
746 std::dynamic_pointer_cast<
747 MultiRegions::DisContField>(m_exp[0]);
748
752 }
753 else
754 {
756 std::dynamic_pointer_cast<MultiRegions::ExpList>(
757 m_exp[0]);
758
759 tmp = MemoryManager<
760 MultiRegions::ExpList>::AllocateSharedPtr(*tmp2);
761 }
762 }
763 }
764 break;
765 case 2:
766 {
767 if (NumHomogeneousDir == 1)
768 {
770 {
771 if (NewField)
772 {
773 bool dealiasing = false;
774
775 tmp = MemoryManager<
776 MultiRegions::ContField3DHomogeneous1D>::
777 AllocateSharedPtr(m_session,
778 m_exp[0]
779 ->GetHomogeneousBasis()
780 ->GetBasisKey(),
781 m_exp[0]->GetHomoLen(),
782 m_useFFT, dealiasing, m_graph,
783 var,
785 }
786 else
787 {
789 tmp2 = std::dynamic_pointer_cast<
790 MultiRegions::ContField3DHomogeneous1D>(
791 m_exp[0]);
792
793 ASSERTL0(tmp2, "Failed to type cast m_exp[0]");
794 tmp = MemoryManager<
795 MultiRegions::ContField3DHomogeneous1D>::
796 AllocateSharedPtr(*tmp2, m_graph, var);
797 }
798 }
800 {
801 if (NewField)
802 {
803 bool dealiasing = false;
804
805 tmp = MemoryManager<
806 MultiRegions::DisContField3DHomogeneous1D>::
807 AllocateSharedPtr(m_session,
808 m_exp[0]
809 ->GetHomogeneousBasis()
810 ->GetBasisKey(),
811 m_exp[0]->GetHomoLen(),
812 m_useFFT, dealiasing, m_graph,
813 var,
815 }
816 else
817 {
819 tmp2 = std::dynamic_pointer_cast<
820 MultiRegions::DisContField3DHomogeneous1D>(
821 m_exp[0]);
822 ASSERTL0(tmp2, "Failed to type cast m_exp[0]");
823
824 tmp = MemoryManager<
825 MultiRegions::DisContField3DHomogeneous1D>::
826 AllocateSharedPtr(*tmp2);
827 }
828 }
829 else
830 {
831 if (NewField)
832 {
833 bool dealiasing = false;
834
835 tmp = MemoryManager<
836 MultiRegions::ExpList3DHomogeneous1D>::
837 AllocateSharedPtr(m_session,
838 m_exp[0]
839 ->GetHomogeneousBasis()
840 ->GetBasisKey(),
841 m_exp[0]->GetHomoLen(),
842 m_useFFT, dealiasing, m_graph,
843 "DefaultVar",
845 }
846 else
847 {
849 std::dynamic_pointer_cast<
850 MultiRegions::ExpList3DHomogeneous1D>(
851 m_exp[0]);
852 ASSERTL0(tmp2, "Failed to type cast m_exp[0]");
853
854 tmp = MemoryManager<
855 MultiRegions::ExpList3DHomogeneous1D>::
856 AllocateSharedPtr(*tmp2);
857 }
858 }
859 }
860 else
861 {
863 {
864 if (NewField)
865 {
868 false,
870 }
871 else // call copy constructor
872 {
873
875 std::dynamic_pointer_cast<
876 MultiRegions::ContField>(m_exp[0]);
877
879 AllocateSharedPtr(*tmp2, m_graph, var);
880 }
881 }
883 {
884 if (NewField)
885 {
888 true,
890 }
891 else // call copy constructor
892 {
894 std::dynamic_pointer_cast<
895 MultiRegions::DisContField>(m_exp[0]);
896
898 AllocateSharedPtr(*tmp2, m_graph, var);
899 }
900 }
901 else
902 {
904 std::dynamic_pointer_cast<MultiRegions::ExpList>(
905 m_exp[0]);
906
907 tmp = MemoryManager<
908 MultiRegions::ExpList>::AllocateSharedPtr(*tmp2);
909 }
910 }
911 }
912 break;
913 case 3:
914 {
916 {
917 if (NewField)
918 {
921 false,
923 }
924 else
925 {
927 std::dynamic_pointer_cast<MultiRegions::ContField>(
928 m_exp[0]);
929
930 tmp = MemoryManager<
931 MultiRegions::ContField>::AllocateSharedPtr(*tmp2,
932 m_graph,
933 var);
934 }
935 }
937 {
938 if (NewField)
939 {
943 }
944 else
945 {
947 std::dynamic_pointer_cast<
948 MultiRegions::DisContField>(m_exp[0]);
949
951 AllocateSharedPtr(*tmp2, m_graph, var);
952 }
953 }
954 else
955 {
957 std::dynamic_pointer_cast<MultiRegions::ExpList>(
958 m_exp[0]);
959
960 tmp =
962 *tmp2);
963 }
964 }
965 break;
966 default:
967 ASSERTL0(false, "Expansion dimension not recognised");
968 break;
969 }
970
971 return tmp;
972 }
#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:341
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:270
SpatialDomains::MeshGraphSharedPtr m_graph
Definition: Field.hpp:102
LibUtilities::SessionReaderSharedPtr m_session
Definition: Field.hpp:101
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 974 of file Field.hpp.

975 {
979 m_exp.clear();
980 m_fielddef = std::vector<LibUtilities::FieldDefinitionsSharedPtr>();
981 m_data = std::vector<std::vector<NekDouble>>();
982 m_variables.clear();
983 }
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:176
std::vector< std::string > m_variables
Definition: Field.hpp:84
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.

◆ CreateExp()

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

Definition at line 125 of file Field.hpp.

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

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

645 {
647 std::string fmt = LibUtilities::FieldIO::GetFileType(filename, c);
648
649 auto it = m_fld.find(fmt);
650
651 if (it == m_fld.end())
652 {
655 m_fld[fmt] = fld;
656 return fld;
657 }
658 else
659 {
660 return it->second;
661 }
662 }
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:327
FieldIOFactory & GetFieldIOFactory()
Returns the FieldIO factory.
Definition: FieldIO.cpp:72
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:57
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:988

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

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

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

326 {
327
329
330 // Set up expansion list
331 int expdim = m_graph->GetMeshDimension();
332 bool dealiasing = false;
333
334 m_session->MatchSolverInfo("USEFFT", "FFTW", m_useFFT, false);
335
336 switch (expdim)
337 {
338 case 1:
339 {
340 ASSERTL0(NumHomogeneousDir <= 2,
341 "Quasi-3D approach is only set up for 1 or 2 "
342 "homogeneous directions");
343
344 if (NumHomogeneousDir == 1)
345 {
347
348 // Define Homogeneous expansion
349 int nplanes;
350 NekDouble ly;
352
353 if (fldfilegiven)
354 {
355 nplanes = m_fielddef[0]->m_numModes[1];
356 ly = m_fielddef[0]->m_homogeneousLengths[0];
357 btype = m_fielddef[0]->m_basis[1];
358 }
359 else
360 {
361 m_session->LoadParameter("HomModesZ", nplanes);
362 m_session->LoadParameter("LY", ly);
364 }
365
366 // Choose points to be at evenly spaced points at
367 // nplanes points
368 const LibUtilities::PointsKey Pkey(
370
371 const LibUtilities::BasisKey Bkey(btype, nplanes, Pkey);
372
375 {
376 ASSERTL0(false, "ContFieldHomogeneous1D or "
377 "DisContFieldHomogenenous1D has "
378 "not been implemented");
379 }
380
381 Exp2DH1 =
384 dealiasing, m_graph,
386 exp = Exp2DH1;
387 }
388 else if (NumHomogeneousDir == 2)
389 {
391
392 int nylines, nzlines;
393 NekDouble ly, lz;
394 LibUtilities::BasisType btype1, btype2;
395
396 if (fldfilegiven)
397 {
398 nylines = m_fielddef[0]->m_numModes[1];
399 nzlines = m_fielddef[0]->m_numModes[2];
400 ly = m_fielddef[0]->m_homogeneousLengths[0];
401 lz = m_fielddef[0]->m_homogeneousLengths[1];
402 btype1 = m_fielddef[0]->m_basis[1];
403 btype2 = m_fielddef[0]->m_basis[2];
404 }
405 else
406 {
407 m_session->LoadParameter("HomModesY", nylines);
408 m_session->LoadParameter("HomModesZ", nzlines);
409 m_session->LoadParameter("LY", ly);
410 m_session->LoadParameter("LZ", lz);
411 btype1 = LibUtilities::eFourier;
412 btype2 = LibUtilities::eFourier;
413 }
414
415 // Choose points to be at evenly spaced points at
416 // nplanes points
417 const LibUtilities::PointsKey PkeyY(
419 const LibUtilities::BasisKey BkeyY(btype1, nylines, PkeyY);
420
421 const LibUtilities::PointsKey PkeyZ(
423 const LibUtilities::BasisKey BkeyZ(btype2, nzlines, PkeyZ);
424
426 {
427 Exp3DH2 = MemoryManager<
428 MultiRegions::ContField3DHomogeneous2D>::
429 AllocateSharedPtr(m_session, BkeyY, BkeyZ, ly, lz,
430 m_useFFT, dealiasing, m_graph,
431 m_session->GetVariable(0),
433 }
435 {
436 Exp3DH2 = MemoryManager<
437 MultiRegions::DisContField3DHomogeneous2D>::
438 AllocateSharedPtr(m_session, BkeyY, BkeyZ, ly, lz,
439 m_useFFT, dealiasing, m_graph,
440 m_session->GetVariable(0),
442 }
443 else
444 {
445 Exp3DH2 = MemoryManager<
446 MultiRegions::ExpList3DHomogeneous2D>::
447 AllocateSharedPtr(m_session, BkeyY, BkeyZ, ly, lz,
448 m_useFFT, dealiasing, m_graph,
450 }
451
452 exp = Exp3DH2;
453 }
454 else
455 {
457
459 {
462 m_session, m_graph, m_session->GetVariable(0),
463 true, false, Collections::eNoCollection);
464 }
466 {
469 m_session->GetVariable(0), true,
471 }
472 else
473 {
476 "DefaultVar",
478 }
479
480 exp = Exp1D;
481 }
482 }
483 break;
484 case 2:
485 {
486 ASSERTL0(NumHomogeneousDir <= 1,
487 "NumHomogeneousDir is only set up for 1");
488
489 if (NumHomogeneousDir == 1)
490 {
492
493 // Define Homogeneous expansion
494 int nplanes;
495 NekDouble lz;
499
500 if (fldfilegiven)
501 {
502 nplanes = m_fielddef[0]->m_numModes[2];
503 lz = m_fielddef[0]->m_homogeneousLengths[0];
504 btype = m_fielddef[0]->m_basis[2];
505
507 {
509 m_fielddef[0]->m_basis[2] =
511 if (nplanes <= 2)
512 {
513 nplanes = 4;
514 }
515 }
516 else if (btype == LibUtilities::eFourierHalfModeRe &&
517 nplanes == 1)
518 {
520 }
521 }
522 else
523 {
524 m_session->LoadParameter("HomModesZ", nplanes);
525 m_session->LoadParameter("LZ", lz);
527 }
528 // Choose points to be at evenly spaced points at
529 // nplanes points
530 const LibUtilities::PointsKey Pkey(nplanes, ptype);
531
532 const LibUtilities::BasisKey Bkey(btype, nplanes, Pkey);
533
535 {
536 Exp3DH1 = MemoryManager<
537 MultiRegions::ContField3DHomogeneous1D>::
538 AllocateSharedPtr(m_session, Bkey, lz, m_useFFT,
539 dealiasing, m_graph,
540 m_session->GetVariable(0),
542 }
543
545 {
546 Exp3DH1 = MemoryManager<
547 MultiRegions::DisContField3DHomogeneous1D>::
548 AllocateSharedPtr(m_session, Bkey, lz, m_useFFT,
549 dealiasing, m_graph,
550 m_session->GetVariable(0),
552 }
553
554 else
555 {
556 Exp3DH1 = MemoryManager<
557 MultiRegions::ExpList3DHomogeneous1D>::
558 AllocateSharedPtr(m_session, Bkey, lz, m_useFFT,
559 dealiasing, m_graph, "DefaultVar",
561 }
562 exp = Exp3DH1;
563 }
564 else
565 {
567
569 {
572 m_session, m_graph, m_session->GetVariable(0),
573 true, false, Collections::eNoCollection);
574 }
575
577 {
580 m_session->GetVariable(0), true,
582 }
583 else
584 {
587 "DefaultVar",
589 }
590
591 exp = Exp2D;
592 }
593 }
594 break;
595 case 3:
596 {
598
600 {
603 m_session->GetVariable(0), true,
605 }
607 {
610 m_session->GetVariable(0), true,
612 }
613 else
614 {
615 Exp3D =
617 m_session, m_graph, true, "DefaultVar",
619 }
620
621 exp = Exp3D;
622 }
623 break;
624 default:
625 ASSERTL0(false, "Expansion dimension not recognised");
626 break;
627 }
628
629 return exp;
630 };
@ 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 107 of file Field.hpp.

◆ m_bndRegionsToWrite

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

Definition at line 106 of file Field.hpp.

◆ m_comm

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

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

Referenced by AppendExpList(), and SetUpFirstExpList().

◆ m_declareExpansionAsDisContField

bool Nektar::FieldUtils::Field::m_declareExpansionAsDisContField

Definition at line 89 of file Field.hpp.

Referenced by AppendExpList(), and SetUpFirstExpList().

◆ m_defComm

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

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

◆ m_fieldPts

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

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

Referenced by FieldIOForFile().

◆ m_graph

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

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

◆ m_nParts

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

Definition at line 98 of file Field.hpp.

◆ m_numHomogeneousDir

int Nektar::FieldUtils::Field::m_numHomogeneousDir

Definition at line 86 of file Field.hpp.

Referenced by CreateExp().

◆ m_partComm

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

Definition at line 97 of file Field.hpp.

◆ m_requireBoundaryExpansion

bool Nektar::FieldUtils::Field::m_requireBoundaryExpansion

Definition at line 91 of file Field.hpp.

Referenced by AppendExpList().

◆ m_session

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

Definition at line 101 of file Field.hpp.

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

◆ m_useFFT

bool Nektar::FieldUtils::Field::m_useFFT

Definition at line 93 of file Field.hpp.

Referenced by AppendExpList(), and SetUpFirstExpList().

◆ m_variables

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

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

◆ m_writeBndFld

bool Nektar::FieldUtils::Field::m_writeBndFld

Definition at line 105 of file Field.hpp.