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  Timer timer;
72 
73  if (m_f->m_verbose)
74  {
75  if(m_f->m_comm->GetRank() == 0)
76  {
77  cout << "ProcessC0Projection: Projecting field into C0 space..."
78  << endl;
79  timer.Start();
80  }
81  }
82 
83  // ensure not using diagonal preconditioner since tends not to converge fo
84  // mass matrix
85  if(m_f->m_graph->GetMeshDimension() == 3)
86  {
87  if(boost::iequals(m_f->m_session->GetSolverInfo("GLOBALSYSSOLN"),
88  "IterativeStaticCond"))
89  {
90  if(boost::iequals(m_f->m_session->GetSolverInfo("PRECONDITIONER"),
91  "Diagonal"))
92  {
93  m_f->m_session->SetSolverInfo("PRECONDITIONER","LowEnergyBlock");
94  }
95  if(boost::iequals(m_f->m_session->GetSolverInfo("PRECONDITIONER"),
96  "FullLinearSpaceWithDiagonal"))
97  {
98  m_f->m_session->SetSolverInfo("PRECONDITIONER","FullLinearSpaceWithLowEnergyBlock");
99  }
100 
101  if(m_f->m_verbose)
102  {
103  if(m_f->m_comm->GetRank() == 0)
104  {
105  cout << "Resetting diagonal precondition to low energy block " << endl;
106  }
107  }
108  }
109  }
110  bool JustPerformLocToGloMap = m_config["localtoglobalmap"].as<bool>();
111  bool HelmSmoother = (boost::iequals(m_config["helmsmoothing"].as<string>(),"Not Set"))? false:true;
112  int nfields = m_f->m_exp.size();
113  Array<OneD, MultiRegions::ExpListSharedPtr> C0ProjectExp(nfields);
114  if(m_config["usexmlbcs"].as<bool>())
115  {
116  for(int i = 0; i < nfields; ++i)
117  {
118  C0ProjectExp[i] = m_f->m_exp[i];
119  }
120  }
121  else
122  {
123  // generate a C0 expansion field with no boundary conditions.
124  bool savedef = m_f->m_declareExpansionAsContField;
125  m_f->m_declareExpansionAsContField = true;
126  C0ProjectExp[0] = m_f->AppendExpList(m_f->m_fielddef[0]->m_numHomogeneousDir,
127  "DefaultVar",true);
128  m_f->m_declareExpansionAsContField = savedef;
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 i =0; i < dim; ++i)
197  {
198  Velocity[i] = 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 
246  if(m_f->m_verbose)
247  {
248  if(m_f->m_comm->GetRank() == 0)
249  {
250  timer.Stop();
251  NekDouble cpuTime = timer.TimePerTest(1);
252 
253  stringstream ss;
254  ss << cpuTime << "s";
255  cout << "C0 Projection CPU Time: " << setw(8) << left
256  << ss.str() << endl;
257  cpuTime = 0.0;
258  }
259  }
260 }
261 }
262 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 Stop()
Definition: Timer.cpp:62
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 Start()
Definition: Timer.cpp:51
NekDouble TimePerTest(unsigned int n)
Returns amount of seconds per iteration in a test with n iterations.
Definition: Timer.cpp:108
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