Nektar++
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// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Computes C0 projection.
32//
33////////////////////////////////////////////////////////////////////////////////
34
35#include <iostream>
36#include <string>
37using namespace std;
38
39#include <boost/core/ignore_unused.hpp>
40
43
44#include "ProcessC0Projection.h"
45
46namespace Nektar
47{
48namespace FieldUtils
49{
50
54 "Computes C0 projection.");
55
57{
58 m_config["fields"] = ConfigOption(false, "All", "Start field to project");
59 m_config["localtoglobalmap"] = ConfigOption(
60 true, "0", "Just perform a local to global mapping and back");
61 m_config["usexmlbcs"] =
62 ConfigOption(true, "0",
63 "Use boundary conditions given in xml file. Requires all "
64 "projected fields to be defined in xml file");
65 m_config["helmsmoothing"] =
66 ConfigOption(false, "Not Set",
67 "Use a Helmholtz smoother to remove high frequency "
68 "components above specified L");
69
70 f->m_declareExpansionAsContField = true;
71}
72
74{
75}
76
77void ProcessC0Projection::v_Process(po::variables_map &vm)
78{
79 m_f->SetUpExp(vm);
80
81 // Skip in case of empty partition
82 if (m_f->m_exp[0]->GetNumElmts() == 0)
83 {
84 return;
85 }
86
87 // ensure not using diagonal preconditioner since tends not to converge for
88 // mass matrix
89 if (m_f->m_graph->GetMeshDimension() == 3)
90 {
91 if (boost::iequals(m_f->m_session->GetSolverInfo("GLOBALSYSSOLN"),
92 "IterativeStaticCond"))
93 {
94 if (boost::iequals(m_f->m_session->GetSolverInfo("PRECONDITIONER"),
95 "Diagonal"))
96 {
97 m_f->m_session->SetSolverInfo("PRECONDITIONER",
98 "LowEnergyBlock");
99 }
100 if (boost::iequals(m_f->m_session->GetSolverInfo("PRECONDITIONER"),
101 "FullLinearSpaceWithDiagonal"))
102 {
103 m_f->m_session->SetSolverInfo(
104 "PRECONDITIONER", "FullLinearSpaceWithLowEnergyBlock");
105 }
106
107 if (m_f->m_verbose)
108 {
109 if (m_f->m_comm->GetRank() == 0)
110 {
111 cout << "Resetting diagonal precondition to low energy "
112 "block "
113 << endl;
114 }
115 }
116 }
117 }
118 bool JustPerformLocToGloMap = m_config["localtoglobalmap"].as<bool>();
119 bool HelmSmoother =
120 (boost::iequals(m_config["helmsmoothing"].as<string>(), "Not Set"))
121 ? false
122 : true;
123 int nfields = m_f->m_exp.size();
125 if (m_config["usexmlbcs"].as<bool>())
126 {
127 for (int i = 0; i < nfields; ++i)
128 {
129 C0ProjectExp[i] = m_f->m_exp[i];
130 }
131 }
132 else
133 {
134 // generate a C0 expansion field with no boundary conditions.
135 bool savedef = m_f->m_declareExpansionAsContField;
136 bool savedef2 = m_f->m_requireBoundaryExpansion;
137 m_f->m_declareExpansionAsContField = true;
138 m_f->m_requireBoundaryExpansion = false;
139 C0ProjectExp[0] =
140 m_f->AppendExpList(m_f->m_numHomogeneousDir, "DefaultVar", true);
141 m_f->m_declareExpansionAsContField = savedef;
142 m_f->m_requireBoundaryExpansion = savedef2;
143 for (int i = 1; i < nfields; ++i)
144 {
145 C0ProjectExp[i] = C0ProjectExp[0];
146 }
147 }
148
149 string fields = m_config["fields"].as<string>();
150 vector<unsigned int> processFields;
152
153 if (fields.compare("All") == 0)
154 {
155 for (int i = 0; i < nfields; ++i)
156 {
157 processFields.push_back(i);
158 }
159 }
160 else
161 {
162 ASSERTL0(ParseUtils::GenerateVector(fields, processFields),
163 "Failed to interpret field string in C0Projection");
164 }
165
166 for (int i = 0; i < processFields.size(); ++i)
167 {
168 ASSERTL0(processFields[i] < nfields,
169 "Attempt to process field that is larger than then number of "
170 "fields available");
171
172 if (m_f->m_verbose)
173 {
174 if (m_f->m_comm->GetRank() == 0)
175 {
176 cout << "\t Processing field: " << processFields[i] << endl;
177 }
178 }
179
180 if (JustPerformLocToGloMap)
181 {
182 int ncoeffs = m_f->m_exp[0]->GetNcoeffs();
183 Array<OneD, NekDouble> Coeffs(ncoeffs);
184 Vmath::Vcopy(ncoeffs, m_f->m_exp[processFields[i]]->GetCoeffs(), 1,
185 Coeffs, 1);
186 C0ProjectExp[processFields[i]]->LocalToGlobal(Coeffs, Coeffs);
187 C0ProjectExp[processFields[i]]->GlobalToLocal(Coeffs, Coeffs);
188
189 Vmath::Vcopy(ncoeffs, Coeffs, 1,
190 tmp = m_f->m_exp[processFields[i]]->UpdateCoeffs(), 1);
191 }
192 else
193 {
194 int ncoeffs = m_f->m_exp[0]->GetNcoeffs();
195 if (HelmSmoother)
196 {
197 int npoints = m_f->m_exp[0]->GetNpoints();
198 NekDouble lambda = m_config["helmsmoothing"].as<NekDouble>();
199 lambda = 2 * M_PI / lambda;
200 lambda = lambda * lambda;
201
202 if (m_f->m_verbose)
203 {
204 cout << "Setting up Helmholtz smoother with lambda = "
205 << lambda << endl;
206 }
207
209 Array<OneD, NekDouble> forcing(npoints);
211
212 Vmath::Smul(npoints, -lambda,
213 m_f->m_exp[processFields[i]]->GetPhys(), 1, forcing,
214 1);
215 Vmath::Zero(ncoeffs,
216 m_f->m_exp[processFields[i]]->UpdateCoeffs(), 1);
217 C0ProjectExp[processFields[i]]->HelmSolve(
218 forcing, m_f->m_exp[processFields[i]]->UpdateCoeffs(),
219 factors);
220 }
221 else
222 {
223 Vmath::Zero(ncoeffs,
224 m_f->m_exp[processFields[i]]->UpdateCoeffs(), 1);
225 C0ProjectExp[processFields[i]]->FwdTrans(
226 m_f->m_exp[processFields[i]]->GetPhys(),
227 m_f->m_exp[processFields[i]]->UpdateCoeffs());
228 }
229 }
230 C0ProjectExp[processFields[i]]->BwdTrans(
231 m_f->m_exp[processFields[i]]->GetCoeffs(),
232 tmp = m_f->m_exp[processFields[i]]->UpdatePhys());
233 }
234}
235} // namespace FieldUtils
236} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
FieldSharedPtr m_f
Field object.
Definition: Module.h:234
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:263
virtual void v_Process(po::variables_map &vm) override
Write mesh to output file.
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
Abstract base class for processing modules.
Definition: Module.h:292
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
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:131
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:991
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:317
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:49
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:408
StdRegions::ConstFactorMap factors
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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:245
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
Represents a command-line configuration option.
Definition: Module.h:131