Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | List of all members
Nektar::FieldUtils::ProcessC0Projection Class Reference

This processing module calculates the Q Criterion and adds it as an extra-field to the output file. More...

#include <ProcessC0Projection.h>

Inheritance diagram for Nektar::FieldUtils::ProcessC0Projection:
[legend]

Public Member Functions

 ProcessC0Projection (FieldSharedPtr f)
 
virtual ~ProcessC0Projection ()
 
virtual void Process (po::variables_map &vm)
 Write mesh to output file. More...
 
virtual std::string GetModuleName ()
 
virtual std::string GetModuleDescription ()
 
virtual ModulePriority GetModulePriority ()
 
- Public Member Functions inherited from Nektar::FieldUtils::ProcessModule
 ProcessModule ()
 
 ProcessModule (FieldSharedPtr p_f)
 
- Public Member Functions inherited from Nektar::FieldUtils::Module
FIELD_UTILS_EXPORT Module (FieldSharedPtr p_f)
 
virtual ~Module ()=default
 
const ConfigOptionGetConfigOption (const std::string &key) const
 
FIELD_UTILS_EXPORT void RegisterConfig (std::string key, std::string value="")
 Register a configuration option with a module. More...
 
FIELD_UTILS_EXPORT void PrintConfig ()
 Print out all configuration options for a module. More...
 
FIELD_UTILS_EXPORT void SetDefaults ()
 Sets default configuration options for those which have not been set. More...
 
FIELD_UTILS_EXPORT void AddFile (std::string fileType, std::string fileName)
 
FIELD_UTILS_EXPORT void EvaluateTriFieldAtEquiSpacedPts (LocalRegions::ExpansionSharedPtr &exp, const Array< OneD, const NekDouble > &infield, Array< OneD, NekDouble > &outfield)
 

Static Public Member Functions

static std::shared_ptr< Modulecreate (FieldSharedPtr f)
 Creates an instance of this class. More...
 

Static Public Attributes

static ModuleKey className
 

Additional Inherited Members

- Public Attributes inherited from Nektar::FieldUtils::Module
FieldSharedPtr m_f
 Field object. More...
 
- Protected Member Functions inherited from Nektar::FieldUtils::Module
 Module ()
 
- Protected Attributes inherited from Nektar::FieldUtils::Module
std::map< std::string, ConfigOptionm_config
 List of configuration values. More...
 
std::set< std::string > m_allowedFiles
 List of allowed file formats. More...
 

Detailed Description

This processing module calculates the Q Criterion and adds it as an extra-field to the output file.

Definition at line 49 of file ProcessC0Projection.h.

Constructor & Destructor Documentation

◆ ProcessC0Projection()

Nektar::FieldUtils::ProcessC0Projection::ProcessC0Projection ( FieldSharedPtr  f)

Definition at line 57 of file ProcessC0Projection.cpp.

57  : ProcessModule(f)
58 {
59  m_config["fields"] = ConfigOption(false, "All", "Start field to project");
60  m_config["localtoglobalmap"] = ConfigOption(
61  true, "0", "Just perform a local to global mapping and back");
62  m_config["usexmlbcs"] = ConfigOption(
63  true, "0", "Use boundary conditions given in xml file. Requires all "
64  "projected fields to be defined in xml file");
65  m_config["helmsmoothing"] = ConfigOption(
66  false, "Not Set", "Use a Helmholtz smoother to remove high frequency "
67  "components above specified L");
68 
69  f->m_declareExpansionAsContField = true;
70 }
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:233

References Nektar::FieldUtils::Module::m_config.

◆ ~ProcessC0Projection()

Nektar::FieldUtils::ProcessC0Projection::~ProcessC0Projection ( )
virtual

Definition at line 72 of file ProcessC0Projection.cpp.

73 {
74 }

Member Function Documentation

◆ create()

static std::shared_ptr<Module> Nektar::FieldUtils::ProcessC0Projection::create ( FieldSharedPtr  f)
inlinestatic

Creates an instance of this class.

Definition at line 53 of file ProcessC0Projection.h.

54  {
56  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ GetModuleDescription()

virtual std::string Nektar::FieldUtils::ProcessC0Projection::GetModuleDescription ( )
inlinevirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 70 of file ProcessC0Projection.h.

71  {
72  return "Projecting field into C0 space";
73  }

◆ GetModuleName()

virtual std::string Nektar::FieldUtils::ProcessC0Projection::GetModuleName ( )
inlinevirtual

Implements Nektar::FieldUtils::Module.

Definition at line 65 of file ProcessC0Projection.h.

66  {
67  return "ProcessC0Projection";
68  }

◆ GetModulePriority()

virtual ModulePriority Nektar::FieldUtils::ProcessC0Projection::GetModulePriority ( )
inlinevirtual

Implements Nektar::FieldUtils::Module.

Definition at line 75 of file ProcessC0Projection.h.

76  {
77  return eModifyExp;
78  }

References Nektar::FieldUtils::eModifyExp.

◆ Process()

void Nektar::FieldUtils::ProcessC0Projection::Process ( po::variables_map &  vm)
virtual

Write mesh to output file.

Implements Nektar::FieldUtils::Module.

Definition at line 76 of file ProcessC0Projection.cpp.

77 {
78  m_f->SetUpExp(vm);
79 
80  // Skip in case of empty partition
81  if (m_f->m_exp[0]->GetNumElmts() == 0)
82  {
83  return;
84  }
85 
86  // ensure not using diagonal preconditioner since tends not to converge fo
87  // mass matrix
88  if (m_f->m_graph->GetMeshDimension() == 3)
89  {
90  if (boost::iequals(m_f->m_session->GetSolverInfo("GLOBALSYSSOLN"),
91  "IterativeStaticCond"))
92  {
93  if (boost::iequals(m_f->m_session->GetSolverInfo("PRECONDITIONER"),
94  "Diagonal"))
95  {
96  m_f->m_session->SetSolverInfo("PRECONDITIONER",
97  "LowEnergyBlock");
98  }
99  if (boost::iequals(m_f->m_session->GetSolverInfo("PRECONDITIONER"),
100  "FullLinearSpaceWithDiagonal"))
101  {
102  m_f->m_session->SetSolverInfo(
103  "PRECONDITIONER", "FullLinearSpaceWithLowEnergyBlock");
104  }
105 
106  if (m_f->m_verbose)
107  {
108  if (m_f->m_comm->GetRank() == 0)
109  {
110  cout << "Resetting diagonal precondition to low energy "
111  "block "
112  << endl;
113  }
114  }
115  }
116  }
117  bool JustPerformLocToGloMap = m_config["localtoglobalmap"].as<bool>();
118  bool HelmSmoother =
119  (boost::iequals(m_config["helmsmoothing"].as<string>(), "Not Set"))
120  ? false
121  : true;
122  int nfields = m_f->m_exp.size();
123  Array<OneD, MultiRegions::ExpListSharedPtr> C0ProjectExp(nfields);
124  if (m_config["usexmlbcs"].as<bool>())
125  {
126  for (int i = 0; i < nfields; ++i)
127  {
128  C0ProjectExp[i] = m_f->m_exp[i];
129  }
130  }
131  else
132  {
133  // generate a C0 expansion field with no boundary conditions.
134  bool savedef = m_f->m_declareExpansionAsContField;
135  bool savedef2 = m_f->m_requireBoundaryExpansion;
136  m_f->m_declareExpansionAsContField = true;
137  m_f->m_requireBoundaryExpansion = false;
138  C0ProjectExp[0] = m_f->AppendExpList(
139  m_f->m_numHomogeneousDir, "DefaultVar", true);
140  m_f->m_declareExpansionAsContField = savedef;
141  m_f->m_requireBoundaryExpansion = savedef2;
142  for (int i = 1; i < nfields; ++i)
143  {
144  C0ProjectExp[i] = C0ProjectExp[0];
145  }
146  }
147 
148  string fields = m_config["fields"].as<string>();
149  vector<unsigned int> processFields;
150 
151  if (fields.compare("All") == 0)
152  {
153  for (int i = 0; i < nfields; ++i)
154  {
155  processFields.push_back(i);
156  }
157  }
158  else
159  {
160  ASSERTL0(
161  ParseUtils::GenerateVector(fields, processFields),
162  "Failed to interpret field string in C0Projection");
163  }
164 
165  for (int i = 0; i < processFields.size(); ++i)
166  {
167  ASSERTL0(processFields[i] < nfields,
168  "Attempt to process field that is larger than then number of "
169  "fields available");
170 
171  if (m_f->m_verbose)
172  {
173  if (m_f->m_comm->GetRank() == 0)
174  {
175  cout << "\t Processing field: " << processFields[i] << endl;
176  }
177  }
178 
179  if (JustPerformLocToGloMap)
180  {
181  int ncoeffs = m_f->m_exp[0]->GetNcoeffs();
182  Vmath::Vcopy(ncoeffs, m_f->m_exp[processFields[i]]->GetCoeffs(), 1,
183  C0ProjectExp[processFields[i]]->UpdateCoeffs(), 1);
184  C0ProjectExp[processFields[i]]->LocalToGlobal();
185  C0ProjectExp[processFields[i]]->GlobalToLocal();
186  Vmath::Vcopy(ncoeffs, C0ProjectExp[processFields[i]]->GetCoeffs(),
187  1, m_f->m_exp[processFields[i]]->UpdateCoeffs(), 1);
188  }
189  else
190  {
191  int ncoeffs = m_f->m_exp[0]->GetNcoeffs();
192  if (HelmSmoother)
193  {
194  int npoints = m_f->m_exp[0]->GetNpoints();
195  NekDouble lambda = m_config["helmsmoothing"].as<NekDouble>();
196  lambda = 2 * M_PI / lambda;
197  lambda = lambda * lambda;
198 
199  if (m_f->m_verbose)
200  {
201  cout << "Setting up Helmholtz smoother with lambda = "
202  << lambda << endl;
203  }
204 
206  Array<OneD, NekDouble> forcing(npoints);
207  factors[StdRegions::eFactorLambda] = lambda;
208 
209  Vmath::Smul(npoints, -lambda,
210  m_f->m_exp[processFields[i]]->GetPhys(), 1, forcing,
211  1);
212 
213  Vmath::Zero(ncoeffs,m_f->m_exp[processFields[i]]->UpdateCoeffs(),1);
214 
215  C0ProjectExp[processFields[i]]->HelmSolve(forcing,
216  m_f->m_exp[processFields[i]]->UpdateCoeffs(), factors);
217  }
218  else
219  {
220  Vmath::Zero(ncoeffs,m_f->m_exp[processFields[i]]->UpdateCoeffs(),1);
221  C0ProjectExp[processFields[i]]->FwdTrans(
222  m_f->m_exp[processFields[i]]->GetPhys(),
223  m_f->m_exp[processFields[i]]->UpdateCoeffs());
224  }
225  }
226  C0ProjectExp[processFields[i]]->BwdTrans(
227  m_f->m_exp[processFields[i]]->GetCoeffs(),
228  m_f->m_exp[processFields[i]]->UpdatePhys());
229  }
230 
231 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
FieldSharedPtr m_f
Field object.
Definition: Module.h:230
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:135
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:314
double NekDouble
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199

References ASSERTL0, Nektar::StdRegions::eFactorLambda, Nektar::ParseUtils::GenerateVector(), Nektar::FieldUtils::Module::m_config, Nektar::FieldUtils::Module::m_f, Vmath::Smul(), Vmath::Vcopy(), and Vmath::Zero().

Member Data Documentation

◆ className

ModuleKey Nektar::FieldUtils::ProcessC0Projection::className
static
Initial value:
=
ModuleKey(eProcessModule, "C0Projection"),
"Computes C0 projection.")
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:290
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:49

Definition at line 57 of file ProcessC0Projection.h.