Nektar++
Functions
FldAddScalGrad.cpp File Reference
#include <cstdio>
#include <cstdlib>
#include <MultiRegions/ExpList.h>
#include <MultiRegions/ExpList1D.h>
#include <MultiRegions/ExpList2D.h>
#include <MultiRegions/ExpList3D.h>
#include <MultiRegions/ContField3D.h>
#include <MultiRegions/ContField2D.h>
Include dependency graph for FldAddScalGrad.cpp:

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Function Documentation

int main ( int  argc,
char *  argv[] 
)

Definition at line 13 of file FldAddScalGrad.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::LibUtilities::SessionReader::CreateInstance(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::LibUtilities::Import(), Nektar::SpatialDomains::MeshGraph::Read(), sign, Vmath::Smul(), Vmath::Svtsvtp(), Vmath::Svtvp(), Vmath::Vvtvp(), Vmath::Vvtvvtp(), and Nektar::LibUtilities::Write().

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