proteus  1.8.1
C/C++/Fortran libraries
Pres.h
Go to the documentation of this file.
1 #ifndef PRES_H
2 #define PRES_H
3 #include <cmath>
4 #include <iostream>
5 #include "CompKernel.h"
6 #include "ModelFactory.h"
7 #include "ArgumentsDict.h"
8 #include "xtensor-python/pyarray.hpp"
9 
10 namespace py = pybind11;
11 
12 namespace proteus
13 {
15  {
16  public:
17  virtual ~cppPres_base(){}
18  virtual void calculateResidual(arguments_dict& args)=0;
19  virtual void calculateJacobian(arguments_dict& args)=0;
20  };
21 
22  template<class CompKernelType,
23  int nSpace,
24  int nQuadraturePoints_element,
25  int nDOF_mesh_trial_element,
26  int nDOF_trial_element,
27  int nDOF_test_element,
28  int nQuadraturePoints_elementBoundary>
29  class cppPres : public cppPres_base
30  {
31  public:
32  CompKernelType ck;
33  cppPres():ck()
34  {}
35  inline
36  void evaluateCoefficients(const double& p,
37  const double& p_last,
38  const double& p_inc,
39  const double massFlux[nSpace],
40  double f[nSpace],
41  double& r)
42  {
43  for (int I=0;I<nSpace;I++)
44  f[I] = massFlux[I];
45  r = p - p_last - p_inc;
46  }
47 
48  inline void calculateElementResidual(//element
49  double* mesh_trial_ref,
50  double* mesh_grad_trial_ref,
51  double* mesh_dof,
52  int* mesh_l2g,
53  double* dV_ref,
54  double* u_trial_ref,
55  double* u_grad_trial_ref,
56  double* u_test_ref,
57  double* u_grad_test_ref,
58  //element boundary
59  double* mesh_trial_trace_ref,
60  double* mesh_grad_trial_trace_ref,
61  double* dS_ref,
62  double* u_trial_trace_ref,
63  double* u_grad_trial_trace_ref,
64  double* u_test_trace_ref,
65  double* u_grad_test_trace_ref,
66  double* normal_ref,
67  double* boundaryJac_ref,
68  //physics
69  int nElements_global,
70  int* u_l2g,
71  double* u_dof,
72  double* q_u,
73  double* q_grad_u,
74  double* q_p_last,
75  double* q_p_inc,
76  double* q_massFlux,
77  double* ebqe_massFlux,
78  double* ebqe_u,
79  double* ebqe_grad_u,
80  int offset_u, int stride_u,
81  double* elementResidual_u,
82  int nExteriorElementBoundaries_global,
83  int* exteriorElementBoundariesArray,
84  int* elementBoundaryElementsArray,
85  int* elementBoundaryLocalElementBoundariesArray,
86  double* element_u,
87  int eN)
88  {
89  for (int i=0;i<nDOF_test_element;i++)
90  {
91  elementResidual_u[i]=0.0;
92  }//i
93  //loop over quadrature points and compute integrands
94  for (int k=0;k<nQuadraturePoints_element;k++)
95  {
96  //compute indeces and declare local storage
97  register int eN_k = eN*nQuadraturePoints_element+k,
98  eN_k_nSpace = eN_k*nSpace;
99  //eN_nDOF_trial_element = eN*nDOF_trial_element;
100  register double u=0.0,grad_u[nSpace],
101  f[nSpace],
102  r=0.0,
103  jac[nSpace*nSpace],
104  jacDet,
105  jacInv[nSpace*nSpace],
106  u_grad_trial[nDOF_trial_element*nSpace],
107  u_test_dV[nDOF_trial_element],
108  u_grad_test_dV[nDOF_test_element*nSpace],
109  dV,x,y,z,
110  G[nSpace*nSpace],G_dd_G,tr_G;
111  //
112  //compute solution and gradients at quadrature points
113  //
114  ck.calculateMapping_element(eN,
115  k,
116  mesh_dof,
117  mesh_l2g,
118  mesh_trial_ref,
119  mesh_grad_trial_ref,
120  jac,
121  jacDet,
122  jacInv,
123  x,y,z);
124  //get the physical integration weight
125  dV = fabs(jacDet)*dV_ref[k];
126  ck.calculateG(jacInv,G,G_dd_G,tr_G);
127  //get the trial function gradients
128  ck.gradTrialFromRef(&u_grad_trial_ref[k*nDOF_trial_element*nSpace],jacInv,u_grad_trial);
129  //get the solution
130  ck.valFromElementDOF(element_u,&u_trial_ref[k*nDOF_trial_element],u);
131  //get the solution gradients
132  ck.gradFromElementDOF(element_u,u_grad_trial,grad_u);
133  //precalculate test function products with integration weights
134  for (int j=0;j<nDOF_trial_element;j++)
135  {
136  u_test_dV[j] = u_test_ref[k*nDOF_trial_element+j]*dV;
137  for (int I=0;I<nSpace;I++)
138  {
139  u_grad_test_dV[j*nSpace+I] = u_grad_trial[j*nSpace+I]*dV;//cek warning won't work for Petrov-Galerkin
140  }
141  }
142  //
143  //calculate pde coefficients at quadrature points
144  //
145  evaluateCoefficients(u,//u is p
146  q_p_last[eN_k],
147  q_p_inc[eN_k],
148  &q_massFlux[eN_k_nSpace],
149  f,
150  r);
151  //
152  //update element residual
153  //
154  for(int i=0;i<nDOF_test_element;i++)
155  {
156  register int i_nSpace=i*nSpace;
157  elementResidual_u[i] += ck.Advection_weak(f,&u_grad_test_dV[i_nSpace]) +
158  ck.Reaction_weak(r,u_test_dV[i]);
159  }//i
160  q_u[eN_k] = u;
161  for(int I=0;I<nSpace;I++)
162  q_grad_u[eN_k_nSpace+I] = grad_u[I];
163  }
164  }
166  {
167  xt::pyarray<double>& mesh_trial_ref = args.array<double>("mesh_trial_ref");
168  xt::pyarray<double>& mesh_grad_trial_ref = args.array<double>("mesh_grad_trial_ref");
169  xt::pyarray<double>& mesh_dof = args.array<double>("mesh_dof");
170  xt::pyarray<int>& mesh_l2g = args.array<int>("mesh_l2g");
171  xt::pyarray<double>& dV_ref = args.array<double>("dV_ref");
172  xt::pyarray<double>& u_trial_ref = args.array<double>("u_trial_ref");
173  xt::pyarray<double>& u_grad_trial_ref = args.array<double>("u_grad_trial_ref");
174  xt::pyarray<double>& u_test_ref = args.array<double>("u_test_ref");
175  xt::pyarray<double>& u_grad_test_ref = args.array<double>("u_grad_test_ref");
176  xt::pyarray<double>& mesh_trial_trace_ref = args.array<double>("mesh_trial_trace_ref");
177  xt::pyarray<double>& mesh_grad_trial_trace_ref = args.array<double>("mesh_grad_trial_trace_ref");
178  xt::pyarray<double>& dS_ref = args.array<double>("dS_ref");
179  xt::pyarray<double>& u_trial_trace_ref = args.array<double>("u_trial_trace_ref");
180  xt::pyarray<double>& u_grad_trial_trace_ref = args.array<double>("u_grad_trial_trace_ref");
181  xt::pyarray<double>& u_test_trace_ref = args.array<double>("u_test_trace_ref");
182  xt::pyarray<double>& u_grad_test_trace_ref = args.array<double>("u_grad_test_trace_ref");
183  xt::pyarray<double>& normal_ref = args.array<double>("normal_ref");
184  xt::pyarray<double>& boundaryJac_ref = args.array<double>("boundaryJac_ref");
185  int nElements_global = args.scalar<int>("nElements_global");
186  xt::pyarray<int>& u_l2g = args.array<int>("u_l2g");
187  xt::pyarray<double>& u_dof = args.array<double>("u_dof");
188  xt::pyarray<double>& q_u = args.array<double>("q_u");
189  xt::pyarray<double>& q_grad_u = args.array<double>("q_grad_u");
190  xt::pyarray<double>& q_p_last = args.array<double>("q_p_last");
191  xt::pyarray<double>& q_p_inc = args.array<double>("q_p_inc");
192  xt::pyarray<double>& q_massFlux = args.array<double>("q_massFlux");
193  xt::pyarray<double>& ebqe_massFlux = args.array<double>("ebqe_massFlux");
194  xt::pyarray<double>& ebqe_u = args.array<double>("ebqe_u");
195  xt::pyarray<double>& ebqe_grad_u = args.array<double>("ebqe_grad_u");
196  int offset_u = args.scalar<int>("offset_u");
197  int stride_u = args.scalar<int>("stride_u");
198  xt::pyarray<double>& globalResidual = args.array<double>("globalResidual");
199  int nExteriorElementBoundaries_global = args.scalar<int>("nExteriorElementBoundaries_global");
200  xt::pyarray<int>& exteriorElementBoundariesArray = args.array<int>("exteriorElementBoundariesArray");
201  xt::pyarray<int>& elementBoundaryElementsArray = args.array<int>("elementBoundaryElementsArray");
202  xt::pyarray<int>& elementBoundaryLocalElementBoundariesArray = args.array<int>("elementBoundaryLocalElementBoundariesArray");
203  //
204  //loop over elements to compute volume integrals and load them into element and global residual
205  //
206  //eN is the element index
207  //eN_k is the quadrature point index for a scalar
208  //eN_k_nSpace is the quadrature point index for a vector
209  //eN_i is the element test function index
210  //eN_j is the element trial function index
211  //eN_k_j is the quadrature point index for a trial function
212  //eN_k_i is the quadrature point index for a trial function
213  for(int eN=0;eN<nElements_global;eN++)
214  {
215  //declare local storage for element residual and initialize
216  register double elementResidual_u[nDOF_test_element],element_u[nDOF_trial_element];
217  for (int i=0;i<nDOF_test_element;i++)
218  {
219  register int eN_i=eN*nDOF_test_element+i;
220  element_u[i] = u_dof.data()[u_l2g.data()[eN_i]];
221  }//i
222  calculateElementResidual(mesh_trial_ref.data(),
223  mesh_grad_trial_ref.data(),
224  mesh_dof.data(),
225  mesh_l2g.data(),
226  dV_ref.data(),
227  u_trial_ref.data(),
228  u_grad_trial_ref.data(),
229  u_test_ref.data(),
230  u_grad_test_ref.data(),
231  mesh_trial_trace_ref.data(),
232  mesh_grad_trial_trace_ref.data(),
233  dS_ref.data(),
234  u_trial_trace_ref.data(),
235  u_grad_trial_trace_ref.data(),
236  u_test_trace_ref.data(),
237  u_grad_test_trace_ref.data(),
238  normal_ref.data(),
239  boundaryJac_ref.data(),
240  nElements_global,
241  u_l2g.data(),
242  u_dof.data(),
243  q_u.data(),
244  q_grad_u.data(),
245  q_p_last.data(),
246  q_p_inc.data(),
247  q_massFlux.data(),
248  ebqe_massFlux.data(),
249  ebqe_u.data(),
250  ebqe_grad_u.data(),
251  offset_u,stride_u,
252  elementResidual_u,
253  nExteriorElementBoundaries_global,
254  exteriorElementBoundariesArray.data(),
255  elementBoundaryElementsArray.data(),
256  elementBoundaryLocalElementBoundariesArray.data(),
257  element_u,
258  eN);
259  //
260  //load element into global residual and save element residual
261  //
262  for(int i=0;i<nDOF_test_element;i++)
263  {
264  register int eN_i=eN*nDOF_test_element+i;
265  globalResidual.data()[offset_u+stride_u*u_l2g.data()[eN_i]]+=elementResidual_u[i];
266  }//i
267  }//elements
268  //
269  //loop over exterior element boundaries to calculate levelset gradient
270  //
271  //ebNE is the Exterior element boundary INdex
272  //ebN is the element boundary INdex
273  //eN is the element index
274  for (int ebNE = 0; ebNE < nExteriorElementBoundaries_global; ebNE++)
275  {
276  register int ebN = exteriorElementBoundariesArray.data()[ebNE],
277  eN = elementBoundaryElementsArray.data()[ebN*2+0],
278  ebN_local = elementBoundaryLocalElementBoundariesArray.data()[ebN*2+0];
279  //eN_nDOF_trial_element = eN*nDOF_trial_element;
280  register double elementResidual_u[nDOF_test_element];
281  double element_u[nDOF_trial_element];
282  for (int i=0;i<nDOF_test_element;i++)
283  {
284  register int eN_i=eN*nDOF_test_element+i;
285  element_u[i] = u_dof.data()[u_l2g.data()[eN_i]];
286  elementResidual_u[i] = 0.0;
287  }//i
288  for (int kb=0;kb<nQuadraturePoints_elementBoundary;kb++)
289  {
290  register int ebNE_kb = ebNE*nQuadraturePoints_elementBoundary+kb,
291  ebNE_kb_nSpace = ebNE_kb*nSpace,
292  ebN_local_kb = ebN_local*nQuadraturePoints_elementBoundary+kb,
293  ebN_local_kb_nSpace = ebN_local_kb*nSpace;
294  register double u_ext=0.0,
295  adv_flux_ext=0.0,
296  grad_u_ext[nSpace],
297  jac_ext[nSpace*nSpace],
298  jacDet_ext,
299  jacInv_ext[nSpace*nSpace],
300  boundaryJac[nSpace*(nSpace-1)],
301  metricTensor[(nSpace-1)*(nSpace-1)],
302  metricTensorDetSqrt,
303  dS,
304  u_test_dS[nDOF_test_element],
305  u_grad_trial_trace[nDOF_trial_element*nSpace],
306  normal[nSpace],x_ext,y_ext,z_ext,
307  G[nSpace*nSpace],G_dd_G,tr_G;
308  //
309  //calculate the solution and gradients at quadrature points
310  //
311  ck.calculateMapping_elementBoundary(eN,
312  ebN_local,
313  kb,
314  ebN_local_kb,
315  mesh_dof.data(),
316  mesh_l2g.data(),
317  mesh_trial_trace_ref.data(),
318  mesh_grad_trial_trace_ref.data(),
319  boundaryJac_ref.data(),
320  jac_ext,
321  jacDet_ext,
322  jacInv_ext,
323  boundaryJac,
324  metricTensor,
325  metricTensorDetSqrt,
326  normal_ref.data(),
327  normal,
328  x_ext,y_ext,z_ext);
329  dS = metricTensorDetSqrt*dS_ref.data()[kb];
330  //get the metric tensor
331  //cek todo use symmetry
332  ck.calculateG(jacInv_ext,G,G_dd_G,tr_G);
333  //compute shape and solution information
334  //shape
335  ck.gradTrialFromRef(&u_grad_trial_trace_ref.data()[ebN_local_kb_nSpace*nDOF_trial_element],jacInv_ext,u_grad_trial_trace);
336  //solution and gradients
337  ck.valFromElementDOF(element_u,&u_trial_trace_ref.data()[ebN_local_kb*nDOF_test_element],u_ext);
338  ck.gradFromElementDOF(element_u,u_grad_trial_trace,grad_u_ext);
339  //precalculate test function products with integration weights
340  for (int j=0;j<nDOF_trial_element;j++)
341  {
342  u_test_dS[j] = u_test_trace_ref.data()[ebN_local_kb*nDOF_test_element+j]*dS;
343  }
344  //
345  //load the boundary values
346  //
347  ebqe_u.data()[ebNE_kb] = u_ext;
348  for (int I=0;I<nSpace;I++)
349  ebqe_grad_u.data()[ebNE_kb_nSpace+I] = grad_u_ext[I];
350  adv_flux_ext=0.0;
351  for (int I=0;I<nSpace;I++)
352  adv_flux_ext += ebqe_massFlux.data()[ebNE_kb_nSpace+I]*normal[I];
353  //
354  //update residuals
355  //
356  for (int i=0;i<nDOF_test_element;i++)
357  {
358  elementResidual_u[i] += ck.ExteriorElementBoundaryFlux(adv_flux_ext,u_test_dS[i]);
359  }//i
360  }//kb
361  //
362  //update the element and global residual storage
363  //
364  for (int i=0;i<nDOF_test_element;i++)
365  {
366  int eN_i = eN*nDOF_test_element+i;
367  globalResidual.data()[offset_u+stride_u*u_l2g.data()[eN_i]] += elementResidual_u[i];
368  }//i
369  }//ebNE
370  }
371 
372  inline void calculateElementJacobian(//element
373  double* mesh_trial_ref,
374  double* mesh_grad_trial_ref,
375  double* mesh_dof,
376  int* mesh_l2g,
377  double* dV_ref,
378  double* u_trial_ref,
379  double* u_grad_trial_ref,
380  double* u_test_ref,
381  double* u_grad_test_ref,
382  //element boundary
383  double* mesh_trial_trace_ref,
384  double* mesh_grad_trial_trace_ref,
385  double* dS_ref,
386  double* u_trial_trace_ref,
387  double* u_grad_trial_trace_ref,
388  double* u_test_trace_ref,
389  double* u_grad_test_trace_ref,
390  double* normal_ref,
391  double* boundaryJac_ref,
392  //physics
393  int nElements_global,
394  double* elementJacobian_u_u,
395  int eN)
396  {
397  for (int i=0;i<nDOF_test_element;i++)
398  for (int j=0;j<nDOF_trial_element;j++)
399  {
400  elementJacobian_u_u[i*nDOF_trial_element+j]=0.0;
401  }
402  for (int k=0;k<nQuadraturePoints_element;k++)
403  {
404  //declare local storage
405  register double jac[nSpace*nSpace],
406  jacDet,
407  jacInv[nSpace*nSpace],
408  dV,
409  u_test_dV[nDOF_test_element],
410  x,y,z;
411  //
412  //calculate solution and gradients at quadrature points
413  //
414  ck.calculateMapping_element(eN,
415  k,
416  mesh_dof,
417  mesh_l2g,
418  mesh_trial_ref,
419  mesh_grad_trial_ref,
420  jac,
421  jacDet,
422  jacInv,
423  x,y,z);
424  //get the physical integration weight
425  dV = fabs(jacDet)*dV_ref[k];
426  for (int j=0;j<nDOF_trial_element;j++)
427  {
428  u_test_dV[j] = u_test_ref[k*nDOF_trial_element+j]*dV;
429  }
430  for(int i=0;i<nDOF_test_element;i++)
431  {
432  for(int j=0;j<nDOF_trial_element;j++)
433  {
434  elementJacobian_u_u[i*nDOF_trial_element+j] += ck.ReactionJacobian_weak(1.0,
435  u_trial_ref[k*nDOF_trial_element+j],
436  u_test_dV[i]);
437  }//j
438  }//i
439  }//k
440  }
442  {
443  xt::pyarray<double>& mesh_trial_ref = args.array<double>("mesh_trial_ref");
444  xt::pyarray<double>& mesh_grad_trial_ref = args.array<double>("mesh_grad_trial_ref");
445  xt::pyarray<double>& mesh_dof = args.array<double>("mesh_dof");
446  xt::pyarray<int>& mesh_l2g = args.array<int>("mesh_l2g");
447  xt::pyarray<double>& dV_ref = args.array<double>("dV_ref");
448  xt::pyarray<double>& u_trial_ref = args.array<double>("u_trial_ref");
449  xt::pyarray<double>& u_grad_trial_ref = args.array<double>("u_grad_trial_ref");
450  xt::pyarray<double>& u_test_ref = args.array<double>("u_test_ref");
451  xt::pyarray<double>& u_grad_test_ref = args.array<double>("u_grad_test_ref");
452  xt::pyarray<double>& mesh_trial_trace_ref = args.array<double>("mesh_trial_trace_ref");
453  xt::pyarray<double>& mesh_grad_trial_trace_ref = args.array<double>("mesh_grad_trial_trace_ref");
454  xt::pyarray<double>& dS_ref = args.array<double>("dS_ref");
455  xt::pyarray<double>& u_trial_trace_ref = args.array<double>("u_trial_trace_ref");
456  xt::pyarray<double>& u_grad_trial_trace_ref = args.array<double>("u_grad_trial_trace_ref");
457  xt::pyarray<double>& u_test_trace_ref = args.array<double>("u_test_trace_ref");
458  xt::pyarray<double>& u_grad_test_trace_ref = args.array<double>("u_grad_test_trace_ref");
459  xt::pyarray<double>& normal_ref = args.array<double>("normal_ref");
460  xt::pyarray<double>& boundaryJac_ref = args.array<double>("boundaryJac_ref");
461  int nElements_global = args.scalar<int>("nElements_global");
462  xt::pyarray<int>& csrRowIndeces_u_u = args.array<int>("csrRowIndeces_u_u");
463  xt::pyarray<int>& csrColumnOffsets_u_u = args.array<int>("csrColumnOffsets_u_u");
464  xt::pyarray<double>& globalJacobian = args.array<double>("globalJacobian");
465  //
466  //loop over elements to compute volume integrals and load them into the element Jacobians and global Jacobian
467  //
468  for(int eN=0;eN<nElements_global;eN++)
469  {
470  register double elementJacobian_u_u[nDOF_test_element*nDOF_trial_element];
471  calculateElementJacobian(mesh_trial_ref.data(),
472  mesh_grad_trial_ref.data(),
473  mesh_dof.data(),
474  mesh_l2g.data(),
475  dV_ref.data(),
476  u_trial_ref.data(),
477  u_grad_trial_ref.data(),
478  u_test_ref.data(),
479  u_grad_test_ref.data(),
480  mesh_trial_trace_ref.data(),
481  mesh_grad_trial_trace_ref.data(),
482  dS_ref.data(),
483  u_trial_trace_ref.data(),
484  u_grad_trial_trace_ref.data(),
485  u_test_trace_ref.data(),
486  u_grad_test_trace_ref.data(),
487  normal_ref.data(),
488  boundaryJac_ref.data(),
489  nElements_global,
490  elementJacobian_u_u,
491  eN);
492  //
493  //load into element Jacobian into global Jacobian
494  //
495  for (int i=0;i<nDOF_test_element;i++)
496  {
497  int eN_i = eN*nDOF_test_element+i;
498  for (int j=0;j<nDOF_trial_element;j++)
499  {
500  int eN_i_j = eN_i*nDOF_trial_element+j;
501  globalJacobian.data()[csrRowIndeces_u_u[eN_i] + csrColumnOffsets_u_u[eN_i_j]] += elementJacobian_u_u[i*nDOF_trial_element+j];
502  }//j
503  }//i
504  }//elements
505  }//computeJacobian
506  };//cppPres
507 
508  inline cppPres_base* newPres(int nSpaceIn,
509  int nQuadraturePoints_elementIn,
510  int nDOF_mesh_trial_elementIn,
511  int nDOF_trial_elementIn,
512  int nDOF_test_elementIn,
513  int nQuadraturePoints_elementBoundaryIn,
514  int CompKernelFlag)
515  {
516  if (nSpaceIn == 2)
517  return proteus::chooseAndAllocateDiscretization2D<cppPres_base,cppPres,CompKernel>(nSpaceIn,
518  nQuadraturePoints_elementIn,
519  nDOF_mesh_trial_elementIn,
520  nDOF_trial_elementIn,
521  nDOF_test_elementIn,
522  nQuadraturePoints_elementBoundaryIn,
523  CompKernelFlag);
524  else
525  return proteus::chooseAndAllocateDiscretization<cppPres_base,cppPres,CompKernel>(nSpaceIn,
526  nQuadraturePoints_elementIn,
527  nDOF_mesh_trial_elementIn,
528  nDOF_trial_elementIn,
529  nDOF_test_elementIn,
530  nQuadraturePoints_elementBoundaryIn,
531  CompKernelFlag);
532  }
533 }//proteus
534 #endif
proteus::cppPres::calculateResidual
void calculateResidual(arguments_dict &args)
Definition: Pres.h:165
proteus::cppPres
Definition: Pres.h:30
proteus::cppPres_base::~cppPres_base
virtual ~cppPres_base()
Definition: Pres.h:17
ModelFactory.h
CompKernel.h
proteus::cppPres_base
Definition: Pres.h:15
proteus::arguments_dict::scalar
T & scalar(const std::string &key)
proteus::arguments_dict::array
xt::pyarray< T > & array(const std::string &key)
proteus::cppPres::cppPres
cppPres()
Definition: Pres.h:33
proteus::cppPres_base::calculateResidual
virtual void calculateResidual(arguments_dict &args)=0
proteus::newPres
cppPres_base * newPres(int nSpaceIn, int nQuadraturePoints_elementIn, int nDOF_mesh_trial_elementIn, int nDOF_trial_elementIn, int nDOF_test_elementIn, int nQuadraturePoints_elementBoundaryIn, int CompKernelFlag)
Definition: Pres.h:508
proteus::cppPres::calculateElementResidual
void calculateElementResidual(double *mesh_trial_ref, double *mesh_grad_trial_ref, double *mesh_dof, int *mesh_l2g, double *dV_ref, double *u_trial_ref, double *u_grad_trial_ref, double *u_test_ref, double *u_grad_test_ref, double *mesh_trial_trace_ref, double *mesh_grad_trial_trace_ref, double *dS_ref, double *u_trial_trace_ref, double *u_grad_trial_trace_ref, double *u_test_trace_ref, double *u_grad_test_trace_ref, double *normal_ref, double *boundaryJac_ref, int nElements_global, int *u_l2g, double *u_dof, double *q_u, double *q_grad_u, double *q_p_last, double *q_p_inc, double *q_massFlux, double *ebqe_massFlux, double *ebqe_u, double *ebqe_grad_u, int offset_u, int stride_u, double *elementResidual_u, int nExteriorElementBoundaries_global, int *exteriorElementBoundariesArray, int *elementBoundaryElementsArray, int *elementBoundaryLocalElementBoundariesArray, double *element_u, int eN)
Definition: Pres.h:48
z
Double * z
Definition: Headers.h:49
u
Double u
Definition: Headers.h:89
proteus::cppPres_base::calculateJacobian
virtual void calculateJacobian(arguments_dict &args)=0
proteus::cppPres::ck
CompKernelType ck
Definition: Pres.h:32
proteus::cppPres::calculateJacobian
void calculateJacobian(arguments_dict &args)
Definition: Pres.h:441
proteus::cppPres::calculateElementJacobian
void calculateElementJacobian(double *mesh_trial_ref, double *mesh_grad_trial_ref, double *mesh_dof, int *mesh_l2g, double *dV_ref, double *u_trial_ref, double *u_grad_trial_ref, double *u_test_ref, double *u_grad_test_ref, double *mesh_trial_trace_ref, double *mesh_grad_trial_trace_ref, double *dS_ref, double *u_trial_trace_ref, double *u_grad_trial_trace_ref, double *u_test_trace_ref, double *u_grad_test_trace_ref, double *normal_ref, double *boundaryJac_ref, int nElements_global, double *elementJacobian_u_u, int eN)
Definition: Pres.h:372
proteus
Definition: ADR.h:17
proteus::f
double f(const double &g, const double &h, const double &hZ)
Definition: SW2DCV.h:58
r
Double r
Definition: Headers.h:83
proteus::arguments_dict
Definition: ArgumentsDict.h:70
proteus::cppPres::evaluateCoefficients
void evaluateCoefficients(const double &p, const double &p_last, const double &p_inc, const double massFlux[nSpace], double f[nSpace], double &r)
Definition: Pres.h:36
ArgumentsDict.h