Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
FldAddScalGrad_elmt.cpp
Go to the documentation of this file.
1 #include <cstdio>
2 #include <cstdlib>
3 
4 #include <MultiRegions/ExpList.h>
10 
11 using namespace std;
12 using namespace Nektar;
13 
14 int main(int argc, char *argv[])
15 {
16  int i,j;
17  int surfID;
18 
19  if(argc != 5)
20  {
21  fprintf(stderr,"Usage: FldAddScalGrad meshfile infld outfld BoundaryID\n");
22  exit(1);
23  }
24 
25  surfID = boost::lexical_cast<int>(argv[argc - 1]);
26 
27  argv[argc -1] = argv[argc - 2];
28 
30  = LibUtilities::SessionReader::CreateInstance(argc, argv);
31 
32  //----------------------------------------------
33  // Read in mesh from input file
34  string meshfile(argv[argc-4]);
35  SpatialDomains::MeshGraphSharedPtr graphShPt = SpatialDomains::MeshGraph::Read(vSession);
36  //----------------------------------------------
37 
38  //----------------------------------------------
39  // Import field file.
40  string fieldfile(argv[argc-3]);
41  vector<LibUtilities::FieldDefinitionsSharedPtr> fielddef;
42  vector<vector<NekDouble> > fielddata;
43  LibUtilities::Import(fieldfile,fielddef,fielddata);
44  //----------------------------------------------
45 
46  //----------------------------------------------
47  // Define Expansion
48  int expdim = graphShPt->GetMeshDimension();
49  int nfields = 1;
50  int addfields = 7;
51  Array<OneD, MultiRegions::ExpListSharedPtr> exp(nfields + addfields);
52  MultiRegions::AssemblyMapCGSharedPtr m_locToGlobalMap;
53 
54  switch(expdim)
55  {
56  case 1:
57  {
58  ASSERTL0(false,"Expansion dimension not recognised");
59  }
60  break;
61  case 2:
62  {
63  ASSERTL0(false,"Expansion dimension not recognised");
64  }
65  break;
66  case 3:
67  {
70  ::AllocateSharedPtr(vSession, graphShPt,
71  vSession->GetVariable(0));
72 
73  m_locToGlobalMap = originalfield->GetLocalToGlobalMap();
74 
75  exp[0] = originalfield;
76  for (i=0; i<expdim; i++)
77  {
79  ::AllocateSharedPtr(vSession, graphShPt,"DefaultVar");
80  }
81 
82  for (i=expdim; i<addfields; i++)
83  {
85  ::AllocateSharedPtr(*originalfield, graphShPt,
86  vSession->GetVariable(0));
87  }
88 
89  }
90  break;
91  default:
92  ASSERTL0(false,"Expansion dimension not recognised");
93  break;
94  }
95  //----------------------------------------------
96 
97  //----------------------------------------------
98  // Copy data from field file
99  for(j = 0; j < nfields+addfields; ++j)
100  {
101  for(int i = 0; i < fielddata.size(); ++i)
102  {
103  exp[j]->ExtractDataToCoeffs(fielddef [i],
104  fielddata[i],
105  fielddef [i]->m_fields[0],
106  exp[j]->UpdateCoeffs());
107  }
108 
109  exp[j]->BwdTrans(exp[j]->GetCoeffs(),exp[j]->UpdatePhys());
110  }
111 
112 
113  //----------------------------------------------
114 
115  //----------------------------------------------
116  int n, cnt, elmtid, nq, offset, nt, boundary, nfq;
117  nt = exp[0]->GetNpoints();
118  Array<OneD, Array<OneD, NekDouble> > grad(expdim);
119  Array<OneD, Array<OneD, NekDouble> > fgrad(expdim);
120  Array<OneD, Array<OneD, NekDouble> > values(addfields-expdim);
121 
122  // Set up mapping from Boundary condition to element details.
125  Array<OneD, int> BoundarytoElmtID;
126  Array<OneD, int> BoundarytoTraceID;
129  Array<OneD, NekDouble> outvalues;
130 
131  exp[0]->GetBoundaryToElmtMap(BoundarytoElmtID,BoundarytoTraceID);
132 
133  //get boundary expansions for each field
134  for (i = 0; i<addfields-expdim; i++)
135  {
136  BndExp[i] = exp[i+1+expdim]->GetBndCondExpansions();
137  }
138 
139 
140  exp[0]->PhysDeriv(exp[0]->GetPhys(),exp[1]->UpdatePhys(),
141  exp[2]->UpdatePhys(),exp[3]->UpdatePhys());
142 
143 
144  // loop over the types of boundary conditions
145  for(cnt = n = 0; n < BndExp[0].num_elements(); ++n)
146  {
147  // identify boundary which the user wanted
148  if(n == surfID)
149  {
150  for(i = 0; i < BndExp[0][n]->GetExpSize(); ++i, cnt++)
151  {
152  // find element and face of this expansion.
153  elmtid = BoundarytoElmtID[cnt];
154  elmt = exp[0]->GetExp(elmtid);
155  nq = elmt->GetTotPoints();
156  offset = exp[0]->GetPhys_Offset(elmtid);
157 
158  // Initialise local arrays for the velocity gradients
159  // size of total number of quadrature points for each element (hence local).
160  for(j = 0; j < expdim; ++j)
161  {
162  grad[j] = Array<OneD, NekDouble>(nq);
163  }
164 
165  if(expdim == 2)
166  {
167  }
168  else
169  {
170  for (j = 0; j< addfields-expdim; j++)
171  {
172  values[j] = BndExp[j][n]->UpdateCoeffs() + BndExp[0][n]->GetCoeff_Offset(i);
173  }
174 
175  // Get face 2D expansion from element expansion
176  bc = boost::dynamic_pointer_cast<StdRegions::StdExpansion2D> (BndExp[0][n]->GetExp(i));
177 
178  // Number of face quadrature points
179  nfq = bc->GetTotPoints();
180 
181  //identify boundary of element
182  boundary = BoundarytoTraceID[cnt];
183 
184  //Extract scalar field
185  U = exp[0]->GetPhys() + offset;
186 
187  //Compute gradients
188  elmt->PhysDeriv(U,grad[0],grad[1],grad[2]);
189 
190 
191  for(j = 0; j < expdim; ++j)
192  {
193  fgrad[j] = Array<OneD, NekDouble>(nfq);
194  }
195 
196 
197  // Get gradient at the quadrature points of the face
198  for(j = 0; j < expdim; ++j)
199  {
200  elmt->GetFacePhysVals(boundary,bc,grad[j],fgrad[j]);
201  //bc->FwdTrans(fgrad[j],values[j]);
202  }
203 
204  const SpatialDomains::GeomFactorsSharedPtr m_metricinfo=bc->GetMetricInfo();
205 
207  = elmt->GetFaceNormal(boundary);
208 
209  Array<OneD, NekDouble> gradnorm(nfq);
210 
211  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
212  {
213  Vmath::Vvtvvtp(nfq,normals[0],1,fgrad[0],1,
214  normals[1],1,fgrad[1],1,gradnorm,1);
215  Vmath::Vvtvp (nfq,normals[2],1,fgrad[2],1,gradnorm,1,gradnorm,1);
216  }
217  else
218  {
219  Vmath::Svtsvtp(nfq,normals[0][0],fgrad[0],1,
220  normals[1][0],fgrad[1],1,gradnorm,1);
221  Vmath::Svtvp(nfq,normals[2][0],fgrad[2],1,gradnorm,1,gradnorm,1);
222  }
223 
224  for(j = 0; j<expdim; j++)
225  {
226  bc->FwdTrans(normals[j],values[j]);
227  }
228 
229  //gradient (grad(u) n)
230  Vmath::Smul(nfq,-1.0,gradnorm,1,gradnorm,1);
231  bc->FwdTrans(gradnorm,values[expdim]);
232  }
233  }
234 
235  }
236  else
237  {
238  cnt += BndExp[0][n]->GetExpSize();
239  }
240  }
241 
242  for(j = 0; j< expdim; j++)
243  {
244  exp[j+1]->FwdTrans(exp[j+1]->GetPhys(), exp[j+1]->UpdateCoeffs());
245  }
246 
247  for(j = 0; j < addfields-expdim; ++j)
248  {
249  int ncoeffs = exp[0]->GetNcoeffs();
250  Array<OneD, NekDouble> output(ncoeffs);
251 
252  output=exp[j+1+expdim]->UpdateCoeffs();
253 
254  int nGlobal=m_locToGlobalMap->GetNumGlobalCoeffs();
255  Array<OneD, NekDouble> outarray(nGlobal,0.0);
256 
257  int bndcnt=0;
258 
259  const Array<OneD,const int>& map = m_locToGlobalMap->GetBndCondCoeffsToGlobalCoeffsMap();
260  NekDouble sign;
261 
262  for(int i = 0; i < BndExp[j].num_elements(); ++i)
263  {
264  if(i==surfID)
265  {
266  const Array<OneD,const NekDouble>& coeffs = BndExp[j][i]->GetCoeffs();
267  for(int k = 0; k < (BndExp[j][i])->GetNcoeffs(); ++k)
268  {
269  sign = m_locToGlobalMap->GetBndCondCoeffsToGlobalCoeffsSign(bndcnt);
270  outarray[map[bndcnt++]] = sign * coeffs[k];
271  }
272  }
273  else
274  {
275  bndcnt += BndExp[j][i]->GetNcoeffs();
276  }
277  }
278  m_locToGlobalMap->GlobalToLocal(outarray,output);
279  }
280 
281 
282  //-----------------------------------------------
283  // Write solution to file with additional computed fields
284  string out(argv[argc-2]);
285  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef
286  = exp[0]->GetFieldDefinitions();
287  std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
288 
289  vector<string > outname;
290 
291  outname.push_back("du/dx");
292  outname.push_back("du/dy");
293  outname.push_back("du/dz");
294  outname.push_back("nx");
295  outname.push_back("ny");
296  outname.push_back("nz");
297  outname.push_back("gradient");
298 
299 
300  for(j = 0; j < nfields+addfields; ++j)
301  {
302  for(i = 0; i < FieldDef.size(); ++i)
303  {
304  if (j >= nfields)
305  {
306  FieldDef[i]->m_fields.push_back(outname[j-nfields]);
307  }
308  else
309  {
310  FieldDef[i]->m_fields.push_back(fielddef[i]->m_fields[j]);
311  }
312  exp[j]->AppendFieldData(FieldDef[i], FieldData[i]);
313  }
314  }
315 
316  LibUtilities::Write(out, FieldDef, FieldData);
317  //-----------------------------------------------
318 
319  return 0;
320 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:22
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:471
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:428
STL namespace.
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
void Import(const std::string &infilename, std::vector< FieldDefinitionsSharedPtr > &fielddefs, std::vector< std::vector< NekDouble > > &fielddata, FieldMetaDataMap &fieldinfomap, const Array< OneD, int > ElementiDs)
Imports an FLD file.
Definition: FieldIO.cpp:115
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
boost::shared_ptr< StdExpansion2D > StdExpansion2DSharedPtr
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< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
Definition: Vmath.cpp:523
void Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
vvtvvtp (scalar times vector plus scalar times vector):
Definition: Vmath.cpp:577
void Write(const std::string &outFile, std::vector< FieldDefinitionsSharedPtr > &fielddefs, std::vector< std::vector< NekDouble > > &fielddata, const FieldMetaDataMap &fieldinfomap)
Write a field file in serial only.
Definition: FieldIO.cpp:81
boost::shared_ptr< ContField3D > ContField3DSharedPtr
Definition: ContField3D.h:191
boost::shared_ptr< AssemblyMapCG > AssemblyMapCGSharedPtr
Definition: AssemblyMapCG.h:52
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:442
Geometry is curved or has non-constant factors.
int main(int argc, char *argv[])