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
41
42#include "ProcessC0Projection.h"
43
44namespace Nektar::FieldUtils
45{
46
50 "Computes C0 projection.");
51
53{
54 m_config["fields"] = ConfigOption(false, "All", "Start field to project");
55 m_config["localtoglobalmap"] = ConfigOption(
56 true, "0", "Just perform a local to global mapping and back");
57 m_config["usexmlbcs"] =
58 ConfigOption(true, "0",
59 "Use boundary conditions given in xml file. Requires all "
60 "projected fields to be defined in xml file");
61 m_config["helmsmoothing"] =
62 ConfigOption(false, "Not Set",
63 "Use a Helmholtz smoother to remove high frequency "
64 "components above specified L");
65
66 f->m_declareExpansionAsContField = true;
67}
68
70{
71}
72
73void ProcessC0Projection::v_Process(po::variables_map &vm)
74{
75 m_f->SetUpExp(vm);
76
77 // Skip in case of empty partition
78 if (m_f->m_exp[0]->GetNumElmts() == 0)
79 {
80 return;
81 }
82
83 // ensure not using diagonal preconditioner since tends not to converge for
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",
94 "LowEnergyBlock");
95 }
96 if (boost::iequals(m_f->m_session->GetSolverInfo("PRECONDITIONER"),
97 "FullLinearSpaceWithDiagonal"))
98 {
99 m_f->m_session->SetSolverInfo(
100 "PRECONDITIONER", "FullLinearSpaceWithLowEnergyBlock");
101 }
102
103 if (m_f->m_verbose)
104 {
105 if (m_f->m_comm->GetRank() == 0)
106 {
107 cout << "Resetting diagonal precondition to low energy "
108 "block "
109 << endl;
110 }
111 }
112 }
113 }
114 bool JustPerformLocToGloMap = m_config["localtoglobalmap"].as<bool>();
115 bool HelmSmoother =
116 (boost::iequals(m_config["helmsmoothing"].as<string>(), "Not Set"))
117 ? false
118 : true;
119 int nfields = m_f->m_exp.size();
121 if (m_config["usexmlbcs"].as<bool>())
122 {
123 for (int i = 0; i < nfields; ++i)
124 {
125 C0ProjectExp[i] = m_f->m_exp[i];
126 }
127 }
128 else
129 {
130 // generate a C0 expansion field with no boundary conditions.
131 bool savedef = m_f->m_declareExpansionAsContField;
132 bool savedef2 = m_f->m_requireBoundaryExpansion;
133 m_f->m_declareExpansionAsContField = true;
134 m_f->m_requireBoundaryExpansion = false;
135 C0ProjectExp[0] =
136 m_f->AppendExpList(m_f->m_numHomogeneousDir, "DefaultVar", true);
137 m_f->m_declareExpansionAsContField = savedef;
138 m_f->m_requireBoundaryExpansion = savedef2;
139 for (int i = 1; i < nfields; ++i)
140 {
141 C0ProjectExp[i] = C0ProjectExp[0];
142 }
143 }
144
145 string fields = m_config["fields"].as<string>();
146 vector<unsigned int> processFields;
148
149 if (fields.compare("All") == 0)
150 {
151 for (int i = 0; i < nfields; ++i)
152 {
153 processFields.push_back(i);
154 }
155 }
156 else
157 {
158 ASSERTL0(ParseUtils::GenerateVector(fields, processFields),
159 "Failed to interpret field string in C0Projection");
160 }
161
162 for (int i = 0; i < processFields.size(); ++i)
163 {
164 ASSERTL0(processFields[i] < nfields,
165 "Attempt to process field that is larger than then number of "
166 "fields available");
167
168 if (m_f->m_verbose)
169 {
170 if (m_f->m_comm->GetRank() == 0)
171 {
172 cout << "\t Processing field: " << processFields[i] << endl;
173 }
174 }
175
176 if (JustPerformLocToGloMap)
177 {
178 int ncoeffs = m_f->m_exp[0]->GetNcoeffs();
179 Array<OneD, NekDouble> Coeffs(ncoeffs);
180 Vmath::Vcopy(ncoeffs, m_f->m_exp[processFields[i]]->GetCoeffs(), 1,
181 Coeffs, 1);
182 C0ProjectExp[processFields[i]]->LocalToGlobal(Coeffs, Coeffs);
183 C0ProjectExp[processFields[i]]->GlobalToLocal(Coeffs, Coeffs);
184
185 Vmath::Vcopy(ncoeffs, Coeffs, 1,
186 tmp = m_f->m_exp[processFields[i]]->UpdateCoeffs(), 1);
187 }
188 else
189 {
190 int ncoeffs = m_f->m_exp[0]->GetNcoeffs();
191 if (HelmSmoother)
192 {
193 int npoints = m_f->m_exp[0]->GetNpoints();
194 NekDouble lambda = m_config["helmsmoothing"].as<NekDouble>();
195 lambda = 2 * M_PI / lambda;
196 lambda = lambda * lambda;
197
198 if (m_f->m_verbose)
199 {
200 cout << "Setting up Helmholtz smoother with lambda = "
201 << lambda << endl;
202 }
203
205 Array<OneD, NekDouble> forcing(npoints);
207
208 Vmath::Smul(npoints, -lambda,
209 m_f->m_exp[processFields[i]]->GetPhys(), 1, forcing,
210 1);
211 Vmath::Zero(ncoeffs,
212 m_f->m_exp[processFields[i]]->UpdateCoeffs(), 1);
213 C0ProjectExp[processFields[i]]->HelmSolve(
214 forcing, m_f->m_exp[processFields[i]]->UpdateCoeffs(),
215 factors);
216 }
217 else
218 {
219 Vmath::Zero(ncoeffs,
220 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 tmp = m_f->m_exp[processFields[i]]->UpdatePhys());
229 }
230}
231} // namespace Nektar::FieldUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
FieldSharedPtr m_f
Field object.
Definition: Module.h:239
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:272
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:301
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
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:130
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:990
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:47
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:402
StdRegions::ConstFactorMap factors
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.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
Represents a command-line configuration option.
Definition: Module.h:129