Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Static Public Member Functions | Static Public Attributes | List of all members
Nektar::Utilities::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::Utilities::ProcessC0Projection:
Inheritance graph
[legend]
Collaboration diagram for Nektar::Utilities::ProcessC0Projection:
Collaboration graph
[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 ()
 
- Public Member Functions inherited from Nektar::Utilities::ProcessModule
 ProcessModule ()
 
 ProcessModule (FieldSharedPtr p_f)
 
 ProcessModule (MeshSharedPtr p_m)
 
- Public Member Functions inherited from Nektar::Utilities::Module
 Module (FieldSharedPtr p_f)
 
void RegisterConfig (string key, string value)
 Register a configuration option with a module. More...
 
void PrintConfig ()
 Print out all configuration options for a module. More...
 
void SetDefaults ()
 Sets default configuration options for those which have not been set. More...
 
bool GetRequireEquiSpaced (void)
 
void SetRequireEquiSpaced (bool pVal)
 
void EvaluateTriFieldAtEquiSpacedPts (LocalRegions::ExpansionSharedPtr &exp, const Array< OneD, const NekDouble > &infield, Array< OneD, NekDouble > &outfield)
 
 Module (MeshSharedPtr p_m)
 
virtual void Process ()=0
 
void RegisterConfig (std::string key, std::string value)
 
void PrintConfig ()
 
void SetDefaults ()
 
MeshSharedPtr GetMesh ()
 
virtual void ProcessVertices ()
 Extract element vertices. More...
 
virtual void ProcessEdges (bool ReprocessEdges=true)
 Extract element edges. More...
 
virtual void ProcessFaces (bool ReprocessFaces=true)
 Extract element faces. More...
 
virtual void ProcessElements ()
 Generate element IDs. More...
 
virtual void ProcessComposites ()
 Generate composites. More...
 
virtual void ClearElementLinks ()
 

Static Public Member Functions

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

Static Public Attributes

static ModuleKey className
 

Additional Inherited Members

- Protected Member Functions inherited from Nektar::Utilities::Module
 Module ()
 
void ReorderPrisms (PerMap &perFaces)
 Reorder node IDs so that prisms and tetrahedra are aligned correctly. More...
 
void PrismLines (int prism, PerMap &perFaces, std::set< int > &prismsDone, std::vector< ElementSharedPtr > &line)
 
- Protected Attributes inherited from Nektar::Utilities::Module
FieldSharedPtr m_f
 Field object. More...
 
map< string, ConfigOptionm_config
 List of configuration values. More...
 
bool m_requireEquiSpaced
 
MeshSharedPtr m_mesh
 Mesh object. More...
 
std::map< std::string,
ConfigOption
m_config
 List of configuration values. More...
 

Detailed Description

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

Definition at line 50 of file ProcessC0Projection.h.

Constructor & Destructor Documentation

Nektar::Utilities::ProcessC0Projection::ProcessC0Projection ( FieldSharedPtr  f)

Definition at line 55 of file ProcessC0Projection.cpp.

References Nektar::Utilities::Module::m_config.

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

Definition at line 65 of file ProcessC0Projection.cpp.

66 {
67 }

Member Function Documentation

static boost::shared_ptr<Module> Nektar::Utilities::ProcessC0Projection::create ( FieldSharedPtr  f)
inlinestatic

Creates an instance of this class.

Definition at line 54 of file ProcessC0Projection.h.

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

55  {
57  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
virtual std::string Nektar::Utilities::ProcessC0Projection::GetModuleName ( )
inlinevirtual

Implements Nektar::Utilities::Module.

Definition at line 66 of file ProcessC0Projection.h.

67  {
68  return "ProcessC0Projection";
69  }
void Nektar::Utilities::ProcessC0Projection::Process ( po::variables_map &  vm)
virtual

Write mesh to output file.

Implements Nektar::Utilities::Module.

Definition at line 69 of file ProcessC0Projection.cpp.

References ASSERTL0, Nektar::StdRegions::eFactorLambda, Nektar::ParseUtils::GenerateOrderedVector(), Nektar::Utilities::Module::m_config, Nektar::Utilities::Module::m_f, Vmath::Smul(), and Vmath::Vcopy().

70 {
71  if (m_f->m_verbose)
72  {
73  if(m_f->m_comm->TreatAsRankZero())
74  {
75  cout << "ProcessC0Projection: Projecting field into C0 space..."
76  << endl;
77  }
78  }
79 
80  // ensure not using diagonal preconditioner since tends not to converge fo
81  // mass matrix
82  if(m_f->m_graph->GetMeshDimension() == 3)
83  {
84  if(boost::iequals(m_f->m_session->GetSolverInfo("GLOBALSYSSOLN"),
85  "IterativeStaticCond"))
86  {
87  if(boost::iequals(m_f->m_session->GetSolverInfo("PRECONDITIONER"),
88  "Diagonal"))
89  {
90  m_f->m_session->SetSolverInfo("PRECONDITIONER","LowEnergyBlock");
91  }
92  if(boost::iequals(m_f->m_session->GetSolverInfo("PRECONDITIONER"),
93  "FullLinearSpaceWithDiagonal"))
94  {
95  m_f->m_session->SetSolverInfo("PRECONDITIONER","FullLinearSpaceWithLowEnergyBlock");
96  }
97 
98  if(m_f->m_verbose)
99  {
100  if(m_f->m_comm->GetRank() == 0)
101  {
102  cout << "Resetting diagonal precondition to low energy block " << endl;
103  }
104  }
105  }
106  }
107  bool JustPerformLocToGloMap = m_config["localtoglobalmap"].as<bool>();
108  bool HelmSmoother = (boost::iequals(m_config["helmsmoothing"].as<string>(),"Not Set"))? false:true;
109  int nfields = m_f->m_exp.size();
110  Array<OneD, MultiRegions::ExpListSharedPtr> C0ProjectExp(nfields);
111  if(m_config["usexmlbcs"].as<bool>())
112  {
113  for(int i = 0; i < nfields; ++i)
114  {
115  C0ProjectExp[i] = m_f->m_exp[i];
116  }
117  }
118  else
119  {
120  // generate a C0 expansion field with no boundary conditions.
121  bool savedef = m_f->m_declareExpansionAsContField;
122  bool savedef2 = m_f->m_requireBoundaryExpansion;
123  m_f->m_declareExpansionAsContField = true;
124  m_f->m_requireBoundaryExpansion = false;
125  C0ProjectExp[0] = m_f->AppendExpList(m_f->m_fielddef[0]->m_numHomogeneousDir,
126  "DefaultVar",true);
127  m_f->m_declareExpansionAsContField = savedef;
128  m_f->m_requireBoundaryExpansion = savedef2;
129  for(int i = 1; i < nfields; ++i)
130  {
131  C0ProjectExp[i] = C0ProjectExp[0];
132  }
133  }
134 
135  string fields = m_config["fields"].as<string>();
136  vector<unsigned int> processFields;
137 
138  if(fields.compare("All") == 0)
139  {
140  for(int i = 0; i < nfields; ++i)
141  {
142  processFields.push_back(i);
143  }
144  }
145  else
146  {
148  processFields),
149  "Failed to interpret field string in C0Projection");
150  }
151 
152  for (int i = 0; i < processFields.size(); ++i)
153  {
154  ASSERTL0(processFields[i] < nfields,
155  "Attempt to process field that is larger than then number of "
156  "fields available");
157 
158  if (m_f->m_verbose)
159  {
160  if(m_f->m_comm->GetRank() == 0)
161  {
162  cout << "\t Processing field: " << processFields[i] << endl;
163  }
164  }
165 
166  if(JustPerformLocToGloMap)
167  {
168  int ncoeffs = m_f->m_exp[0]->GetNcoeffs();
169  Vmath::Vcopy(ncoeffs,m_f->m_exp[processFields[i]]->GetCoeffs(),1,
170  C0ProjectExp[processFields[i]]->UpdateCoeffs(),1);
171  C0ProjectExp[processFields[i]]->LocalToGlobal();
172  C0ProjectExp[processFields[i]]->GlobalToLocal();
173  Vmath::Vcopy(ncoeffs,C0ProjectExp[processFields[i]]->GetCoeffs(),1,
174  m_f->m_exp[processFields[i]]->UpdateCoeffs(),1);
175  }
176  else
177  {
178  if(HelmSmoother)
179  {
180  int dim = m_f->m_graph->GetSpaceDimension();
181  int npoints = m_f->m_exp[0]->GetNpoints();
182  NekDouble lambda = m_config["helmsmoothing"].as<NekDouble>();
183  lambda = 2*M_PI/lambda;
184  lambda = lambda*lambda;
185 
186  if(m_f->m_verbose)
187  {
188  cout << "Setting up Helmholtz smoother with lambda = " << lambda << endl;
189  }
190 
192  Array<OneD, NekDouble> forcing(npoints);
193  factors[StdRegions::eFactorLambda] = -lambda;
194 
195  Array<OneD, Array<OneD, NekDouble> > Velocity(dim);
196  for(int j =0; j < dim; ++j)
197  {
198  Velocity[j] = Array<OneD, NekDouble> (npoints,0.0);
199  }
200 
201  C0ProjectExp[processFields[i]]->BwdTrans(m_f->m_exp[processFields[i]]->GetCoeffs(),
202  m_f->m_exp[processFields[i]]->UpdatePhys());
203 
204  Vmath::Smul(npoints,-lambda,m_f->m_exp[processFields[i]]->GetPhys(),1,
205  forcing,1);
206 
207  // Note we are using the
208  // LinearAdvectionDiffusionReaction solver here
209  // instead of HelmSolve since lambda is negative and
210  // so matrices are not positive definite. Ideally
211  // should allow for negative lambda coefficient in
212  // HelmSolve
213  C0ProjectExp[processFields[i]]->LinearAdvectionDiffusionReactionSolve(Velocity,
214  forcing,
215  m_f->m_exp[processFields[i]]->UpdateCoeffs(),
216  -lambda);
217  }
218  else
219  {
220  C0ProjectExp[processFields[i]]->BwdTrans(m_f->m_exp[processFields[i]]->GetCoeffs(),
221  m_f->m_exp[processFields[i]]->UpdatePhys());
222  C0ProjectExp[processFields[i]]->FwdTrans(m_f->m_exp[processFields[i]]->GetPhys(),
223  m_f->m_exp[processFields[i]]->UpdateCoeffs());
224  }
225  }
226  }
227 
228  // reset FieldDef in case of serial input and parallel output
229  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef
230  = m_f->m_exp[0]->GetFieldDefinitions();
231  // reset up FieldData with new values before projecting
232  std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
233 
234  for(int i = 0; i < nfields; ++i)
235  {
236  for (int j = 0; j < FieldDef.size(); ++j)
237  {
238  FieldDef[j]->m_fields.push_back(m_f->m_fielddef[0]->m_fields[i]);
239  m_f->m_exp[i]->AppendFieldData(FieldDef[j], FieldData[j]);
240  }
241  }
242 
243  m_f->m_fielddef = FieldDef;
244  m_f->m_data = FieldData;
245 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
static bool GenerateOrderedVector(const char *const str, std::vector< unsigned int > &vec)
Definition: ParseUtils.hpp:97
map< string, ConfigOption > m_config
List of configuration values.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:251
FieldSharedPtr m_f
Field object.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047

Member Data Documentation

ModuleKey Nektar::Utilities::ProcessC0Projection::className
static
Initial value:

Definition at line 58 of file ProcessC0Projection.h.