Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
ProcessC0Projection.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: ProcessC0Projection.cpp
4 //
5 // For more information, please see: http://www.nektar.info/
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Computes C0 projection.
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
36 #include <string>
37 #include <iostream>
38 using namespace std;
39 
40 #include "ProcessC0Projection.h"
41 
44 
45 namespace Nektar
46 {
47 namespace Utilities
48 {
49 
50 ModuleKey ProcessC0Projection::className =
52  ModuleKey(eProcessModule, "C0Projection"),
53  ProcessC0Projection::create, "Computes C0 projection.");
54 
55 ProcessC0Projection::ProcessC0Projection(FieldSharedPtr f) : 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 }
64 
66 {
67 }
68 
69 void ProcessC0Projection::Process(po::variables_map &vm)
70 {
71  if (m_f->m_verbose)
72  {
73  if(m_f->m_comm->GetRank() == 0)
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  m_f->m_declareExpansionAsContField = true;
123  C0ProjectExp[0] = m_f->AppendExpList(m_f->m_fielddef[0]->m_numHomogeneousDir,
124  "DefaultVar",true);
125  m_f->m_declareExpansionAsContField = savedef;
126  for(int i = 1; i < nfields; ++i)
127  {
128  C0ProjectExp[i] = C0ProjectExp[0];
129  }
130  }
131 
132  string fields = m_config["fields"].as<string>();
133  vector<unsigned int> processFields;
134 
135  if(fields.compare("All") == 0)
136  {
137  for(int i = 0; i < nfields; ++i)
138  {
139  processFields.push_back(i);
140  }
141  }
142  else
143  {
145  processFields),
146  "Failed to interpret field string in C0Projection");
147  }
148 
149  for (int i = 0; i < processFields.size(); ++i)
150  {
151  ASSERTL0(processFields[i] < nfields,
152  "Attempt to process field that is larger than then number of "
153  "fields available");
154 
155  if (m_f->m_verbose)
156  {
157  if(m_f->m_comm->GetRank() == 0)
158  {
159  cout << "\t Processing field: " << processFields[i] << endl;
160  }
161  }
162 
163  if(JustPerformLocToGloMap)
164  {
165  int ncoeffs = m_f->m_exp[0]->GetNcoeffs();
166  Vmath::Vcopy(ncoeffs,m_f->m_exp[processFields[i]]->GetCoeffs(),1,
167  C0ProjectExp[processFields[i]]->UpdateCoeffs(),1);
168  C0ProjectExp[processFields[i]]->LocalToGlobal();
169  C0ProjectExp[processFields[i]]->GlobalToLocal();
170  Vmath::Vcopy(ncoeffs,C0ProjectExp[processFields[i]]->GetCoeffs(),1,
171  m_f->m_exp[processFields[i]]->UpdateCoeffs(),1);
172  }
173  else
174  {
175  if(HelmSmoother)
176  {
177  int dim = m_f->m_graph->GetSpaceDimension();
178  int npoints = m_f->m_exp[0]->GetNpoints();
179  NekDouble lambda = m_config["helmsmoothing"].as<NekDouble>();
180  lambda = 2*M_PI/lambda;
181  lambda = lambda*lambda;
182 
183  if(m_f->m_verbose)
184  {
185  cout << "Setting up Helmholtz smoother with lambda = " << lambda << endl;
186  }
187 
189  Array<OneD, NekDouble> forcing(npoints);
190  factors[StdRegions::eFactorLambda] = -lambda;
191 
192  Array<OneD, Array<OneD, NekDouble> > Velocity(dim);
193  for(int j =0; j < dim; ++j)
194  {
195  Velocity[j] = Array<OneD, NekDouble> (npoints,0.0);
196  }
197 
198  C0ProjectExp[processFields[i]]->BwdTrans(m_f->m_exp[processFields[i]]->GetCoeffs(),
199  m_f->m_exp[processFields[i]]->UpdatePhys());
200 
201  Vmath::Smul(npoints,-lambda,m_f->m_exp[processFields[i]]->GetPhys(),1,
202  forcing,1);
203 
204  // Note we are using the
205  // LinearAdvectionDiffusionReaction solver here
206  // instead of HelmSolve since lambda is negative and
207  // so matrices are not positive definite. Ideally
208  // should allow for negative lambda coefficient in
209  // HelmSolve
210  C0ProjectExp[processFields[i]]->LinearAdvectionDiffusionReactionSolve(Velocity,
211  forcing,
212  m_f->m_exp[processFields[i]]->UpdateCoeffs(),
213  -lambda);
214  }
215  else
216  {
217  C0ProjectExp[processFields[i]]->BwdTrans(m_f->m_exp[processFields[i]]->GetCoeffs(),
218  m_f->m_exp[processFields[i]]->UpdatePhys());
219  C0ProjectExp[processFields[i]]->FwdTrans(m_f->m_exp[processFields[i]]->GetPhys(),
220  m_f->m_exp[processFields[i]]->UpdateCoeffs());
221  }
222  }
223  }
224 
225  // reset FieldDef in case of serial input and parallel output
226  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef
227  = m_f->m_exp[0]->GetFieldDefinitions();
228  // reset up FieldData with new values before projecting
229  std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
230 
231  for(int i = 0; i < nfields; ++i)
232  {
233  for (int j = 0; j < FieldDef.size(); ++j)
234  {
235  FieldDef[j]->m_fields.push_back(m_f->m_fielddef[0]->m_fields[i]);
236  m_f->m_exp[i]->AppendFieldData(FieldDef[j], FieldData[j]);
237  }
238  }
239 
240  m_f->m_fielddef = FieldDef;
241  m_f->m_data = FieldData;
242 }
243 }
244 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
pair< ModuleType, string > ModuleKey
static bool GenerateOrderedVector(const char *const str, std::vector< unsigned int > &vec)
Definition: ParseUtils.hpp:97
virtual void Process()=0
map< string, ConfigOption > m_config
List of configuration values.
STL namespace.
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
boost::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:695
Represents a command-line configuration option.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
ModuleFactory & GetModuleFactory()
Abstract base class for processing modules.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215