proteus  1.8.1
C/C++/Fortran libraries
ErrorResidualMethod.cpp
Go to the documentation of this file.
1 #include "MeshAdaptPUMI.h"
2 #include <PCU.h>
3 #include <petscksp.h>
4 
5 #include <apf.h>
6 #include <apfMesh.h>
7 #include <apfShape.h>
8 #include <apfDynamicMatrix.h>
9 #include <apfNumbering.h>
10 
11 #include <iostream>
12 #include <fstream>
13 
19 //Global variables used to make it easier to pass these variables from MeshAdaptPUMIDrvr
20 int approx_order; //shape function order
21 int int_order; //integration order
23 double a_kl = 0.5; //flux term weight
24 
25 inline void getProps(double*rho,double*nu)
26 //Function used to transfer MeshAdaptPUMIDrvr variables into global variables
27 {
28  rho_0 = rho[0];
29  nu_0 = nu[0];
30  rho_1 = rho[1];
31  nu_1 = nu[1];
32  return;
33 }
34 
35 double MeshAdaptPUMIDrvr::getMPvalue(double field_val,double val_0, double val_1)
46 {
47  return val_0*(1-field_val)+val_1*field_val;
48 }
49 
50 apf::Vector3 getFaceNormal(apf::Mesh* mesh, apf::MeshEntity* face)
51 //Function used to get the unit vector normal to an element face
52 {
53  apf::Vector3 normal;
54  apf::Adjacent verts;
55  mesh->getAdjacent(face,0,verts);
56  //apf::Vector3 vtxs[verts.getSize()];
57  apf::Vector3 vtxs[3];
58  for(int i=0;i<verts.getSize();i++){
59  mesh->getPoint(verts[i],0,vtxs[i]);
60  }
61  apf::Vector3 a,b;
62  if(mesh->getDimension()==2){
63  vtxs[2] = vtxs[0];
64  vtxs[2][2] = 1.0;
65  }
66  a = vtxs[1]-vtxs[0];
67  b = vtxs[2]-vtxs[0];
68  normal = apf::cross(a,b);
69 
70  return normal.normalize();
71 }
72 
73 double getDotProduct(apf::Vector3 a, apf::Vector3 b)
74 //Function used to get the dot product between two vectors
75 {
76  return (a[0]*b[0] + a[1]*b[1] + a[2]*b[2]);
77 }
78 
79 double getDotProduct(apf::Matrix3x3 a, apf::Matrix3x3 b)
80 //Overloaded function used to get the dot product analog between two matrices
81 {
82  double temp =0;
83  for(int i=0;i<3;i++){
84  for(int j=0;j<3;j++){
85  temp = temp + a[i][j]*b[i][j];
86  }
87  }
88  return temp;
89 }
90 
91 
92 bool isInTet(apf::Mesh* mesh, apf::MeshEntity* ent, apf::Vector3 pt)
93 //Function used to test if a point pt is inside the tetrahedron ent
94 //Returns a boolean: 1 if is in tet, 0 if not
95 {
96  bool isin=0;
97 
98  apf::Adjacent verts;
99  mesh->getAdjacent(ent,0,verts);
100  apf::Vector3 vtxs[4]; //4 points
101  for(int i=0;i<4;i++){
102  mesh->getPoint(verts[i],0,vtxs[i]);
103  }
104  apf::Vector3 c[4];
105  c[0] = vtxs[1]-vtxs[0];
106  c[1] = vtxs[2]-vtxs[0];
107  c[2] = vtxs[3]-vtxs[0];
108  c[3] = pt-vtxs[0];
109 
110  apf::Matrix3x3 K,Kinv;
111  apf::Vector3 F;
112  for(int i=0;i<3;i++){
113  for(int j=0;j<3;j++){
114  K[i][j] = getDotProduct(c[i],c[j]);
115  }
116  F[i] = getDotProduct(c[3],c[i]);
117  }
118  Kinv = invert(K);
119  apf::DynamicMatrix Kinv_dyn = apf::fromMatrix(Kinv);
120  apf::DynamicVector F_dyn = apf::fromVector(F);
121  apf::DynamicVector uvw; //result
122  apf::multiply(Kinv_dyn,F_dyn,uvw);
123  if(uvw[0] >= 0 && uvw[1] >=0 && uvw[2] >=0 && (uvw[0]+uvw[1]+uvw[2])<=1) isin = 1;
124  return isin;
125 }
126 
127 bool isInSimplex(apf::Mesh* mesh, apf::MeshEntity* ent, apf::Vector3 pt, int dim)
128 //Function used to test if a point pt is inside the tetrahedron ent
129 //Returns a boolean: 1 if is in tet, 0 if not
130 {
131  bool isin=0;
132  int numverts = dim+1;
133  apf::Adjacent verts;
134  mesh->getAdjacent(ent,0,verts);
135  apf::Vector3 vtxs[4];
136  for(int i=0;i<numverts;i++){
137  mesh->getPoint(verts[i],0,vtxs[i]);
138  }
139  apf::Vector3 c[4];
140  if(dim==2){
141  c[0] = vtxs[1]-vtxs[0];
142  c[1] = vtxs[2]-vtxs[0];
143  c[2] = pt-vtxs[0];
144  }
145  else if(dim==3){
146  c[0] = vtxs[1]-vtxs[0];
147  c[1] = vtxs[2]-vtxs[0];
148  c[2] = vtxs[3]-vtxs[0];
149  c[3] = pt-vtxs[0];
150  }
151 
152  apf::Matrix3x3 K(0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0);
153  apf::Matrix3x3 Kinv;
154  apf::Vector3 F(0.0,0.0,0.0);
155  for(int i=0;i<dim;i++){
156  for(int j=0;j<dim;j++){
157  K[i][j] = getDotProduct(c[i],c[j]);
158  }
159  F[i] = getDotProduct(c[dim],c[i]);
160  }
161  if(dim==2)
162  K[2][2] = 1.0;
163  Kinv = invert(K);
164  apf::DynamicMatrix Kinv_dyn = apf::fromMatrix(Kinv);
165  apf::DynamicVector F_dyn = apf::fromVector(F);
166  apf::DynamicVector uvw; //result
167  apf::multiply(Kinv_dyn,F_dyn,uvw);
168  if(dim==2){
169  if(uvw[0] >= 0 && uvw[1] >=0 && (uvw[0]+uvw[1])<=1) isin = 1;
170  }
171  else{
172  if(uvw[0] >= 0 && uvw[1] >=0 && uvw[2] >=0 && (uvw[0]+uvw[1]+uvw[2])<=1) isin = 1;
173  }
174  return isin;
175 }
176 
177 
178 /*
179 double a_k(apf::Matrix3x3 u, apf::Matrix3x3 v,double nu){
180  //u and v are gradients of a vector
181  apf::Matrix3x3 temp_u = u+apf::transpose(u);
182  apf::Matrix3x3 temp_v = v+apf::transpose(v);
183  return nu*getDotProduct(temp_u,temp_v);
184 }
185 
186 double b_k(double a, apf::Matrix3x3 b){
187  //b is a gradient of a vector
188  return a*(b[0][0]+b[1][1]+b[1][1]);
189 }
190 
191 double c_k(apf::Vector3 a, apf::Matrix3x3 b, apf::Vector3 c){
192  //b is a gradient of a vector
193  return getDotProduct(b*a,c);
194 }
195 */
196 
197 void getLHS(Mat &K,apf::NewArray <apf::DynamicVector> &shdrv,int nsd,double weight, double visc_val,int nshl)
198 //Function used to get the LHS of the local error problem.
199 //The LHS is entirely A(\phi,\phi) which can be decomposed into a diagonal contributions and off-diagonal contributions
200 //Inputs:
201 // shdrv is the set of shape function derivatives for an element evaluated at a quadrature point
202 // nsd is the number of spatial dimensions
203 // weight is the corresponding weight for a given quadrature point
204 // visc_val is the viscosity at that quadrature point
205 // nshl is the number of local shape functions in an element
206 //Outputs:
207 // K is the matrix representing the LHS
208 {
209  PetscScalar term1[nshl][nshl], term2[nshl][nshl];
210  //Calculate LHS Diagonal Block Term
211  for(int s=0; s<nshl;s++){
212  for(int t=0; t<nshl;t++){
213  double temp=0;
214  for(int j=0;j<nsd;j++){
215  temp+=shdrv[s][j]*shdrv[t][j];
216  }
217  term1[s][t] = temp*weight*visc_val;
218  }
219  }
220  int idx[nshl]; //indices for PETSc Mat insertion
221  for(int i=0; i<nsd;i++){
222  for(int j=0;j<nshl;j++){
223  idx[j] = i*nshl+j;
224  }
225  MatSetValues(K,nshl,idx,nshl,idx,term1[0],ADD_VALUES);
226  }
227  int idxr[nshl],idxc[nshl]; //indices for PETSc rows and columns
228  for(int i = 0; i< nsd;i++){
229  for(int j=0; j< nsd;j++){
230  for(int s=0;s<nshl;s++){
231  for(int t=0;t<nshl;t++){
232  term2[s][t] = shdrv[s][j]*shdrv[t][i]*weight*visc_val;
233  }
234  }
235  for(int count=0;count<nshl;count++){
236  idxr[count] = i*nshl+count;
237  idxc[count] = j*nshl+count;
238  }
239  MatSetValues(K,nshl,idxr,nshl,idxc,term2[0],ADD_VALUES);
240  }
241  } //end 2nd term loop
242 }
243 
244 void getRHS(Vec &F,apf::NewArray <double> &shpval,apf::NewArray <apf::DynamicVector> &shdrv,apf::Vector3 vel_vect,apf::Matrix3x3 grad_vel,
245  int nsd,double weight,int nshl,
246  double visc_val,double density,apf::Vector3 grad_density,double pressure,
247  double g[3])
248 //Function used to get the RHS of the local error problem.
249 //The RHS is the weak residual of the N-S equations
250 //Inputs:
251 // shpval are the local shape functions evaluated at a quadrature point
252 // shdrv are the local shape function derivatives evaluated at a quadrature point
253 // vel_vect is the velocity vector at a quadrature point
254 // grad_vel is the velocity gradient at a quadrature point
255 // nsd is the number of spatial dimensions
256 // weight is the corresponding weight for a given quadrature point
257 // nshl is the number of local shape functions in an element
258 // visc_val is the viscosity at a quadrature point
259 // density is the density at a quadrature point
260 // grad_density is the density gradient at a quadrature point
261 // pressure is the pressure at a quadrature point
262 // g is the gravity vector
263 //Outputs:
264 // F is the vector representing the RHS
265 
266 {
267  int idx[nshl];
268  for( int i = 0; i<nsd; i++){
269  double temp_vect[nshl];
270  for( int s=0;s<nshl;s++){
271  idx[s] = i*nshl+s;
272 
273  //forcing term
274  //temp_vect[s] = (g[i]+0.0)*shpval[s];
275  //temp_vect[s] += pressure/density*shdrv[s][i]; //pressure term
276  double force = (g[i]+0.0)*shpval[s];
277  double pressure_force = pressure/density*shdrv[s][i];
278  double a_rho_term = 0;
279  double b_rho_term = -pressure/(density*density)*grad_density[i]*shpval[s];
280  double a_term = 0;
281  double c_term = 0;
282  //a(u,v) and c(u,u,v) term
283  for(int j=0;j<nsd;j++){
284  a_term += -visc_val*shdrv[s][j]*(grad_vel[i][j]+grad_vel[j][i]);
285  a_rho_term += visc_val*(grad_vel[i][j]+grad_vel[j][i])*shpval[s]*grad_density[j]/(density);
286  c_term += -shpval[s]*grad_vel[i][j]*vel_vect[j];
287  }
288  temp_vect[s] = force+pressure_force+a_term+c_term;
289  //temp_vect[s] = force+pressure_force+a_rho_term+b_rho_term+a_term+c_term;
290  temp_vect[s] = temp_vect[s]*weight;
291  } //end loop over number of shape functions
292  VecSetValues(F,nshl,idx,temp_vect,ADD_VALUES);
293  } //end loop over spatial dimensions
294 }
295 
296 
297 void MeshAdaptPUMIDrvr::computeDiffusiveFlux(apf::Mesh*m,apf::Field* voff, apf::Field* visc,apf::Field* pref, apf::Field* velf)
312 {
313  if(comm_rank==0)
314  std::cerr<<"Begin computeDiffusiveFlux()"<<std::endl;
315  int numbqpt, nshl;
316  int hier_off;
317  if(nsd==2)
318  hier_off=3;
319  else if(nsd==3)
320  hier_off=4;
321  apf::MeshEntity* bent,*ent;
322  apf::MeshIterator* iter = m->begin(nsd-1); //loop over faces
323 
324  //Need to get number of bqpt
325  while(bent = m->iterate(iter)){
326  apf::MeshElement* b_elem;
327  b_elem = apf::createMeshElement(m,bent);
328  numbqpt = apf::countIntPoints(b_elem,int_order);
329  apf::destroyMeshElement(b_elem);
330  break;
331  }
332  m->end(iter);
333 
334  diffFlux = m->createDoubleTag("diffFlux",numbqpt*nsd*2);
335  apf::MeshElement* tempelem; apf::Element * tempvelo,*temppres,*tempvoff;
336  apf::MeshElement* b_elem;
337  apf::Adjacent adjFaces;
338  apf::Vector3 normal,centerdir;
339  int orientation;
340  double tempflux[numbqpt*nsd];
341  double *flux; flux = (double*) calloc(numbqpt*nsd*2,sizeof(double));
342  apf::NewArray <double> shpval;
343  apf::NewArray <double> shpval_temp;
344 
345  apf::FieldShape* err_shape = apf::getHierarchic(2);
346  apf::EntityShape* elem_shape;
347  apf::Vector3 bqpt,bqptl,bqptshp;
348  double weight, Jdet;
349  apf::Matrix3x3 J;
350  apf::Matrix3x3 tempgrad_velo;
351  apf::Matrix3x3 identity(1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0);
352 
353  iter=m->begin(nsd-1);
354  while(ent=m->iterate(iter))
355  {
356  m->setDoubleTag(ent,diffFlux,flux);
357  }
358  m->end(iter);
359  free(flux);
360 
361  if(comm_rank==0)
362  std::cerr<<"Initialized flux"<<std::endl;
363  //loop over regions
364  PCU_Comm_Begin();
365  iter = m->begin(nsd);
366  int ent_count=0;
367  while(ent = m->iterate(iter))
368  {
369  //Shape functions of the region and not the boundaries
370  if(ent_count==0){
371  nshl=apf::countElementNodes(err_shape,m->getType(ent));
372  shpval_temp.allocate(nshl);
373  nshl= nshl-hier_off;
374  shpval.allocate(nshl);
375  elem_shape = err_shape->getEntityShape(m->getType(ent));
376  }
377 
378  m->getAdjacent(ent,nsd-1,adjFaces);
379  for(int adjcount =0;adjcount<adjFaces.getSize();adjcount++){
380  bent = adjFaces[adjcount];
381  normal=getFaceNormal(m,bent);
382  centerdir=apf::getLinearCentroid(m,ent)-apf::getLinearCentroid(m,bent);
383  //if(isInTet(m,ent,apf::project(normal,centerdir)*centerdir.getLength()+apf::getLinearCentroid(m,bent)))
384  if(isInSimplex(m,ent,apf::project(normal,centerdir).normalize()*centerdir.getLength()+apf::getLinearCentroid(m,bent),nsd)){
385  orientation = 1;
386  }
387  else{
388  orientation = 0;
389  }
390 
391  //begin calculation of flux
392  b_elem = apf::createMeshElement(m,bent);
393  tempelem = apf::createMeshElement(m,ent);
394  temppres = apf::createElement(pref,tempelem);
395  tempvelo = apf::createElement(velf,tempelem);
396  tempvoff = apf::createElement(voff,tempelem);
397 
398  for(int l = 0;l<numbqpt;l++)
399  {
400  apf::Vector3 bflux(0.0,0.0,0.0);
401  apf::Matrix3x3 tempbflux(0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0);
402  apf::getIntPoint(b_elem,int_order,l,bqpt);
403  weight = apf::getIntWeight(b_elem,int_order,l);
404  apf::getJacobian(b_elem,bqpt,J); //evaluate the Jacobian at the quadrature point
405  Jdet=fabs(apf::getJacobianDeterminant(J,nsd-1));
406  bqptl=apf::boundaryToElementXi(m,bent,ent,bqpt);
407  apf::getVectorGrad(tempvelo,bqptl,tempgrad_velo);
408  tempgrad_velo = apf::transpose(tempgrad_velo);
409 
410  apf::ModelEntity* me=m->toModel(bent);
411  int tag = m->getModelTag(me);
412  apf::ModelEntity* boundary_face = m->findModelEntity(nsd-1,tag);
413 
414  if(me==boundary_face && has_gBC){
415  int BCtype[4];
416  double fluxdata[4][numbqpt];
417  //for(int i=1;i<nsd+1;i++){ //ignores 0th index because that's pressure
418  // m->getIntTag(bent,BCtag[i],&(BCtype[i]));
419  //}
420  m->getIntTag(bent,BCtag,&(BCtype[0]));
421  if((BCtype[1]+BCtype[2]+BCtype[3] != 3) && BCtype[1] == 1 ){
422  std::cerr << "diffusive flux not fully specified on face " << localNumber(bent) << '\n';
423  std::cerr << "BCtype "<<BCtype[1]<<" "<<BCtype[2]<<" "<<BCtype[3]<<std::endl;
424  abort();
425  }
426  if(BCtype[1]+BCtype[2]+BCtype[3] == 3){
427  for(int i=1;i<nsd+1;i++)
428  m->getDoubleTag(bent,fluxtag[i],&(fluxdata[i][0]));
429  bflux = apf::Vector3(fluxdata[1][l],fluxdata[2][l],fluxdata[3][l]);
430  bflux = bflux-identity*apf::getScalar(temppres,bqptl)/getMPvalue(apf::getScalar(tempvoff,bqptl),rho_0,rho_1)*normal;
431  }
432  else{
433  tempbflux = (tempgrad_velo+apf::transpose(tempgrad_velo))*getMPvalue(apf::getScalar(tempvoff,bqptl),nu_0,nu_1)
434  -identity*apf::getScalar(temppres,bqptl)/getMPvalue(apf::getScalar(tempvoff,bqptl),rho_0,rho_1);
435  bflux = tempbflux*normal;
436  }
437  }
438  else{
439  tempbflux = (tempgrad_velo+apf::transpose(tempgrad_velo))*getMPvalue(apf::getScalar(tempvoff,bqptl),nu_0,nu_1)
440  -identity*apf::getScalar(temppres,bqptl)/getMPvalue(apf::getScalar(tempvoff,bqptl),rho_0,rho_1);
441  bflux = tempbflux*normal;
442  } //end if boundary
443  bflux = bflux*weight*Jdet;
444  bflux.toArray(&(tempflux[l*nsd]));
445  }
446  flux = (double*) calloc(numbqpt*nsd*2,sizeof(double));
447  m->getDoubleTag(bent,diffFlux,flux);
448  for (int i=0;i<numbqpt*nsd;i++){
449  flux[orientation*numbqpt*nsd+i] = tempflux[i];
450  }
451  m->setDoubleTag(bent,diffFlux,flux);
452  free(flux);
453  apf::destroyMeshElement(tempelem);apf::destroyElement(tempvelo);apf::destroyElement(temppres); apf::destroyElement(tempvoff);
454 
455  //Parallel Communications
456  apf::ModelEntity* me=m->toModel(bent);
457  apf::ModelEntity* boundary_face = m->findModelEntity(nsd-1,m->getModelTag(me));
458  apf::Copies remotes;
459  if(m->isShared(bent))
460  {
461  m->getRemotes(bent,remotes);
462  for(apf::Copies::iterator it=remotes.begin(); it!=remotes.end();++it)
463  {
464  PCU_COMM_PACK(it->first, it->second);
465  PCU_COMM_PACK(it->first, orientation);
466  PCU_COMM_PACK(it->first, tempflux);
467  }
468  } //end if
469  ent_count++;
470  } //end loop over faces
471  } //end loop over regions
472  m->end(iter);
473  if(comm_rank==0)
474  std::cerr<<"Sending flux"<<std::endl;
475  PCU_Comm_Send();
476  flux = (double*) calloc(numbqpt*nsd*2,sizeof(double));
477  while(PCU_Comm_Receive())
478  {
479  PCU_COMM_UNPACK(bent);
480  PCU_COMM_UNPACK(orientation);
481  PCU_COMM_UNPACK(tempflux);
482  m->getDoubleTag(bent,diffFlux,flux);
483  for (int i=0;i<numbqpt*nsd;i++){
484  flux[orientation*numbqpt*nsd+i] = flux[orientation*numbqpt*nsd+i]+tempflux[i];
485  }
486  m->setDoubleTag(bent,diffFlux,flux);
487  }
488  PCU_Barrier();
489  free(flux);
490  if(comm_rank==0)
491  std::cerr<<"End computeDiffusiveFlux()"<<std::endl;
492 }
493 
494 void MeshAdaptPUMIDrvr::getBoundaryFlux(apf::Mesh* m, apf::MeshEntity* ent, double * endflux)
502 {
503  int nshl;
504  apf::NewArray <double> shpval;
505  apf::NewArray <double> shpval_temp;
506  double* flux;
507 
508  apf::FieldShape* err_shape = apf::getHierarchic(2);
509  apf::EntityShape* elem_shape;
510 
511  //loop over element faces
512  apf::Adjacent boundaries;
513  apf::MeshEntity* bent;
514  apf::MeshElement* b_elem;
515  apf::Vector3 bqpt,bqptl,bqptshp;
516 
517  apf::Vector3 normal;
518  apf::Vector3 centerdir;
519 
520  //Shape functions of the region and not the boundaries
521  nshl=apf::countElementNodes(err_shape,m->getType(ent));
522  shpval_temp.allocate(nshl);
523  int hier_off;
524  if(nsd==2)
525  hier_off=3;
526  else if(nsd==3)
527  hier_off=4;
528  nshl= nshl-hier_off;
529  shpval.allocate(nshl);
530  elem_shape = err_shape->getEntityShape(m->getType(ent));
531 
532  m->getAdjacent(ent,nsd-1,boundaries);
533  for(int adjcount =0;adjcount<boundaries.getSize();adjcount++){
534 
535  apf::Vector3 bflux(0.0,0.0,0.0);
536  bent = boundaries[adjcount];
537 
538  b_elem = apf::createMeshElement(m,bent);
539  normal=getFaceNormal(m,bent);
540  centerdir=apf::getLinearCentroid(m,ent)-apf::getLinearCentroid(m,bent);
541  int orientation = 0;
542 
543  if(isInSimplex(m,ent,apf::project(normal,centerdir).normalize()*centerdir.getLength()+apf::getLinearCentroid(m,bent),nsd)){
544  normal = normal*-1.0; //normal needs to face the other direction
545  orientation=1;
546  }
547  apf::ModelEntity* me=m->toModel(bent);
548  int tag = m->getModelTag(me);
549  apf::ModelEntity* boundary_face = m->findModelEntity(nsd-1,tag);
550 
551  double flux_weight[2];
552  if(me==boundary_face){
553  if(orientation==0){flux_weight[0]=1; flux_weight[1]=0;}
554  else{flux_weight[0]=0; flux_weight[1]=-1;}
555  }
556  else{
557  if(orientation==0){flux_weight[0]=(1-a_kl); flux_weight[1]=a_kl;}
558  else{ flux_weight[0]=-a_kl; flux_weight[1]=-1*(1-a_kl);}
559  }
560 
561  int numbqpt = apf::countIntPoints(b_elem,int_order);
562  flux = (double*) calloc(numbqpt*nsd*2,sizeof(double));
563  m->getDoubleTag(bent,diffFlux,flux);
564  for(int l=0; l<numbqpt;l++){
565  apf::getIntPoint(b_elem,int_order,l,bqpt);
566  bqptshp=apf::boundaryToElementXi(m,bent,ent,bqpt);
567  elem_shape->getValues(NULL,NULL,bqptshp,shpval_temp);
568  for(int j=0;j<nshl;j++){shpval[j] = shpval_temp[hier_off+j];}
569  for(int i=0;i<nsd;i++){
570  for(int s=0;s<nshl;s++){
571  endflux[i*nshl+s] = endflux[i*nshl+s]+(flux_weight[0]*flux[l*nsd+i]+flux_weight[1]*flux[numbqpt*nsd+l*nsd+i])*shpval[s];
572  }
573  }
574  }//end of boundary integration loop
575  free(flux);
576  } //end for adjacent faces
577 }
578 
579 apf::Field* MeshAdaptPUMIDrvr::getViscosityField(apf::Field* voff)
586 {
587  apf::Field* visc = apf::createLagrangeField(m,"viscosity",apf::SCALAR,1);
588  apf::MeshEntity* ent;
589  apf::MeshIterator* iter = m->begin(0);
590  double vof_val, visc_val;
591  while(ent = m->iterate(iter)){ //loop through all vertices
592  vof_val=apf::getScalar(voff,ent,0);
593  visc_val = getMPvalue(vof_val,nu_0, nu_1);
594  apf::setScalar(visc, ent, 0,visc_val);
595  }
596  m->end(iter);
597  return visc;
598 }
599 
600 
601 void setErrorField(apf::Field* estimate,Vec coef,apf::MeshEntity* ent,int nsd,int nshl)
602 //Function used to store the computed coefficients from the local error problem onto a field
603 {
604  apf::Mesh* m = apf::getMesh(estimate);
605  double coef_ez[nshl*nsd];
606  int ez_idx[nshl*nsd];
607  for(int ez=0;ez<nshl*nsd;ez++){ez_idx[ez]=ez;}
608  VecGetValues(coef,nshl*nsd,ez_idx,coef_ez);
609 
610  //Copy coefficients onto field
611  apf::Adjacent adjvert;
612  m->getAdjacent(ent,0,adjvert);
613  for(int idx=0;idx<adjvert.getSize();idx++){
614  double coef_sub[3]={0,0,0};
615  apf::setVector(estimate,adjvert[idx],0,&coef_sub[0]);
616  }
617 
618  apf::Adjacent adjedg;
619  m->getAdjacent(ent,1,adjedg);
620  for(int idx=0;idx<nshl;idx++){
621  double coef_sub[3] ={coef_ez[idx],coef_ez[nshl+idx],coef_ez[nshl*2+idx]};
622  apf::setVector(estimate,adjedg[idx],0,&coef_sub[0]);
623  }
624 }
625 
633 {
634  if(comm_rank==0) std::cout<<"Start removing BC tags/data"<<std::endl;
635  apf::MeshEntity* ent;
636  apf::MeshIterator* fIter = m->begin(nsd-1);
637  while(ent=m->iterate(fIter))
638  {
639  if(has_gBC && m->hasTag(ent,BCtag)){
640  m->removeTag(ent,BCtag);
641  for(int i=0;i<4;i++)
642  {
643  if(i>0 && m->hasTag(ent,fluxtag[i]))
644  m->removeTag(ent,fluxtag[i]);
645  }
646  }
647  if(m->hasTag(ent,diffFlux))
648  m->removeTag(ent,diffFlux);
649 
650  }
651  m->end(fIter);
652  if(has_gBC){
653  m->destroyTag(BCtag);
654  for(int i=0;i<4;i++)
655  {
656  if(i>0)
657  m->destroyTag(fluxtag[i]);
658  }
659  }
660  m->destroyTag(diffFlux);
661  if(comm_rank==0) std::cerr<<"Destroyed BC and flux tags"<<std::endl;
662 }
663 
664 void MeshAdaptPUMIDrvr::get_local_error(double &total_error)
665 
673 {
674  getProps(rho,nu);
677  nsd = m->getDimension();
678 
679  //***** Get Solution Fields First *****//
680  apf::Field* voff = m->findField("vof");
681  assert(voff);
682  apf::Field* velf = m->findField("velocity");
683  assert(velf);
684  apf::Field* pref = m->findField("p");
685  assert(pref);
686  //***** *****//
687 
688  //***** Compute the viscosity field *****//
689  apf::Field* visc = getViscosityField(voff);
690 
691  //***** Compute diffusive flux *****//
692  computeDiffusiveFlux(m,voff,visc,pref,velf);
693 
694  //Initialize the Error Fields
695  freeField(err_reg);
696  freeField(errRho_reg);
697  freeField(errRel_reg);
698  //err_reg = apf::createField(m,"ErrorRegion",apf::VECTOR,apf::getConstant(nsd));
699  err_reg = apf::createField(m,"ErrorRegion",apf::SCALAR,apf::getVoronoiShape(nsd,1));
700  errRho_reg = apf::createField(m,"ErrorDensity",apf::SCALAR,apf::getVoronoiShape(nsd,1));
701  errRel_reg = apf::createField(m,"RelativeError",apf::SCALAR,apf::getVoronoiShape(nsd,1));
702 
703  //Start computing element quantities
704  int numqpt; //number of quadrature points
705  int nshl; //number of local shape functions
706  int elem_type; //what type of topology
707  double weight; //value container for the weight at each qpt
708  double Jdet;
709  apf::FieldShape* err_shape = apf::getHierarchic(approx_order);
710  apf::Field * estimate = apf::createField(m, "err_est", apf::VECTOR, err_shape);
711  apf::EntityShape* elem_shape;
712  apf::Vector3 qpt; //container for quadrature points
713  apf::MeshElement* element;
714  apf::Element* visc_elem, *pres_elem,*velo_elem,*vof_elem;
715  apf::Element* est_elem;
716  apf::Matrix3x3 J; //actual Jacobian matrix
717  apf::Matrix3x3 invJ; //inverse of Jacobian
718  apf::NewArray <double> shpval; //array to store shape function values at quadrature points
719  apf::NewArray <double> shpval_temp; //array to store shape function values at quadrature points temporarily
720  apf::NewArray <apf::Vector3> shgval; //array to store shape function values at quadrature points
721 
722  apf::DynamicMatrix invJ_copy;
723  apf::NewArray <apf::DynamicVector> shdrv;
724  apf::NewArray <apf::DynamicVector> shgval_copy;
725 
726  apf::MeshIterator* iter = m->begin(nsd); //loop over elements
727  apf::MeshEntity* ent;
728 
729 
730  double err_est = 0;
731  double err_est_total=0;
732  double u_norm_total=0;
733  errRho_max = 0;
734  while(ent = m->iterate(iter)){ //loop through all elements
735 
736  elem_type = m->getType(ent);
737  if(!(elem_type != 4 || elem_type != 2)){ //2|TRI, 4|TET
738  std::cout<<"Not a Tri or Tet present"<<std::endl;
739  exit(0);
740  }
741  element = apf::createMeshElement(m,ent);
742  pres_elem = apf::createElement(pref,element);
743  velo_elem = apf::createElement(velf,element);
744  visc_elem = apf::createElement(visc,element); //at vof currently
745  vof_elem = apf::createElement(voff,element);
746 
747  numqpt=apf::countIntPoints(element,int_order); //generally p*p maximum for shape functions
748  nshl=apf::countElementNodes(err_shape,elem_type);
749  shgval.allocate(nshl);
750  shpval_temp.allocate(nshl);
751 
752  int hier_off;//there is an offset that needs to be made to isolate the hierarchic edge modes
753  if(nsd ==2)
754  hier_off=3;
755  else if(nsd==3)
756  hier_off= 4;
757  nshl = nshl - hier_off;
758 
759  shpval.allocate(nshl); shgval_copy.allocate(nshl); shdrv.allocate(nshl);
760 
761  //LHS Matrix Initialization
762  int ndofs = nshl*nsd;
763  Mat K; //matrix size depends on nshl, which may vary from element to element
764  MatCreate(PETSC_COMM_SELF,&K);
765  MatSetSizes(K,ndofs,ndofs,ndofs,ndofs);
766  MatSetFromOptions(K);
767  MatSetUp(K); //is this inefficient? check later
768 
769  //RHS Vector Initialization
770  Vec F;
771  VecCreate(PETSC_COMM_SELF,&F);
772  VecSetSizes(F,ndofs,ndofs);
773  VecSetUp(F);
774 
775  //loop through all qpts
776  for(int k=0;k<numqpt;k++){
777  apf::getIntPoint(element,int_order,k,qpt); //get a quadrature point and store in qpt
778  apf::getJacobian(element,qpt,J); //evaluate the Jacobian at the quadrature point
779  J = apf::transpose(J); //Is PUMI still defined in this way?
780  if(nsd==2)
781  J[2][2] = 1.0; //this is necessary to avoid singular matrix
782  invJ = invert(J);
783  Jdet=fabs(apf::getJacobianDeterminant(J,nsd));
784  weight = apf::getIntWeight(element,int_order,k);
785  invJ_copy = apf::fromMatrix(invJ);
786 
787  //first get the shape function values for error shape functions
788  elem_shape = err_shape->getEntityShape(elem_type);
789  elem_shape->getValues(NULL,NULL,qpt,shpval_temp);
790  elem_shape->getLocalGradients(NULL,NULL,qpt,shgval);
791 
792  for(int i =0;i<nshl;i++){ //get the true derivative and copy only the edge modes for use
793  shgval_copy[i] = apf::fromVector(shgval[i+hier_off]);
794  shpval[i] = shpval_temp[i+hier_off];
795  apf::multiply(shgval_copy[i],invJ_copy,shdrv[i]);
796  }
797 
798  //obtain needed values
799 
800  apf::Vector3 vel_vect;
801  apf::Matrix3x3 grad_vel;
802  apf::getVector(velo_elem,qpt,vel_vect);
803  apf::getVectorGrad(velo_elem,qpt,grad_vel);
804  grad_vel = apf::transpose(grad_vel);
805  apf::Vector3 grad_vof;
806  apf::getGrad(vof_elem,qpt,grad_vof);
807 
808  double density = getMPvalue(apf::getScalar(vof_elem,qpt),rho_0,rho_1);
809  double pressure = apf::getScalar(pres_elem,qpt);
810  double visc_val = apf::getScalar(visc_elem,qpt);
811  apf::Vector3 grad_rho = grad_vof*(rho_1-rho_0);
812 
813  //Left-Hand Side
814  getLHS(K,shdrv,nsd,weight,visc_val,nshl);
815 
816  //Get RHS
817  getRHS(F,shpval,shdrv,vel_vect,grad_vel,nsd,weight,nshl,visc_val,density,grad_rho,pressure,g);
818 
819  } // end quadrature loop
820 
821  //to complete integration, scale by the determinant of the Jacobian
822 
823  MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
824  MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
825  MatScale(K,Jdet); //must be done after assembly
826  VecAssemblyBegin(F);
827  VecAssemblyEnd(F);
828  VecScale(F,Jdet); //must be done after assembly
829  double* bflux;
830  int F_idx[ndofs];
831  bflux = (double*) calloc(ndofs,sizeof(double));
832 
833  getBoundaryFlux(m, ent,bflux);
834  for(int s=0;s<ndofs;s++){
835  F_idx[s]=s;
836  }
837  VecSetValues(F,ndofs,F_idx,bflux,ADD_VALUES);
838  VecAssemblyBegin(F); VecAssemblyEnd(F);
839  free(bflux);
840  Vec coef;
841  VecCreate(PETSC_COMM_SELF,&coef);
842  VecSetSizes(coef,ndofs,ndofs);
843  VecSetUp(coef);
844 
845  KSP ksp; //initialize solver context
846  KSPCreate(PETSC_COMM_SELF,&ksp);
847  KSPSetOperators(ksp,K,K);
848  KSPSetType(ksp,KSPPREONLY);
849  PC pc;
850  KSPGetPC(ksp,&pc);
851  PCSetType(pc,PCLU);
852  KSPSetFromOptions(ksp);
853 
854  KSPSolve(ksp,F,coef);
855 
856  KSPDestroy(&ksp); //destroy ksp
857 
858  setErrorField(estimate,coef,ent,nsd,nshl);
859 
860  //compute the local error
861  double Acomp=0;
862  double Bcomp=0;
863  double visc_avg=0;
864  double u_norm = 0;
865  apf::Matrix3x3 phi_ij;
866  apf::Matrix3x3 vel_ij;
867  apf::Vector3 vel_vect;
868  apf::Vector3 grad_vof;
869 
870  est_elem= apf::createElement(estimate,element);
871  for(int k=0; k<numqpt;k++){
872  apf::getIntPoint(element,int_order,k,qpt); //get a quadrature point and store in qpt
873  apf::getJacobian(element,qpt,J); //evaluate the Jacobian at the quadrature point
874 
875  invJ = invert(J);
876  invJ = apf::transpose(invJ);
877  Jdet=fabs(apf::getJacobianDeterminant(J,nsd));
878  weight = apf::getIntWeight(element,int_order,k);
879  invJ_copy = apf::fromMatrix(invJ);
880 
881  //first get the shape function values for error shape functions
882  elem_shape = err_shape->getEntityShape(elem_type);
883  elem_shape->getValues(NULL,NULL,qpt,shpval_temp);
884  elem_shape->getLocalGradients(NULL,NULL,qpt,shgval);
885 
886  for(int i =0;i<nshl;i++){ //get the true derivative and copy only the edge modes for use
887  shgval_copy[i] = apf::fromVector(shgval[i+hier_off]);
888  shpval[i] = shpval_temp[i+hier_off];
889  apf::multiply(shgval_copy[i],invJ_copy,shdrv[i]);
890  }
891  double visc_val = apf::getScalar(visc_elem,qpt);
892  double pres_val = apf::getScalar(pres_elem,qpt);
893  double density = getMPvalue(apf::getScalar(vof_elem,qpt),rho_0,rho_1);
894  apf::getVectorGrad(est_elem,qpt,phi_ij);
895  apf::getGrad(vof_elem,qpt,grad_vof);
896  apf::getVector(velo_elem,qpt,vel_vect);
897  apf::getVectorGrad(velo_elem,qpt,vel_ij);
898  vel_ij = apf::transpose(vel_ij);
899  phi_ij = apf::transpose(phi_ij);
900 
901  Acomp = Acomp + visc_val*getDotProduct(phi_ij,phi_ij+apf::transpose(phi_ij))*weight;
902  Bcomp = Bcomp + apf::getDiv(velo_elem,qpt)*apf::getDiv(velo_elem,qpt)*weight;
903  visc_avg = visc_avg + visc_val*weight;
904  u_norm = u_norm + visc_val*getDotProduct(vel_ij,vel_ij+apf::transpose(vel_ij))*weight;
905 
906  } //end compute local error
907  visc_avg = visc_avg*Jdet/apf::measure(element);
908  Acomp = Acomp*Jdet/visc_avg; //nondimensionalize with average viscosity, Jacobians can cancel out, but this is done for clarity
909  Bcomp = Bcomp*Jdet;
910  u_norm = u_norm/visc_avg*Jdet;
911  err_est = sqrt(Acomp);
912 
913  apf::Vector3 err_in(err_est,Acomp,Bcomp);
914  //apf::setVector(err_reg,ent,0,err_in);
915  apf::setScalar(err_reg,ent,0,err_est);
916  double errRho = err_est/sqrt(apf::measure(element));
917  apf::setScalar(errRho_reg,ent,0,errRho);
918  if(errRho>errRho_max)
919  errRho_max = errRho;
920 
921  double err_rel = err_est/sqrt(u_norm);
922  apf::setScalar(errRel_reg,ent,0,err_rel);
923 
924  err_est_total = err_est_total+(Acomp); //for tracking the upper bound
925  u_norm_total = u_norm_total + u_norm;
926 
927  MatDestroy(&K); //destroy the matrix
928  VecDestroy(&F); //destroy vector
929  VecDestroy(&coef); //destroy vector
930 
931  apf::destroyElement(visc_elem);apf::destroyElement(pres_elem);apf::destroyElement(velo_elem);apf::destroyElement(est_elem);apf::destroyElement(vof_elem);
932  } //end element loop
933 
934  PCU_Add_Doubles(&err_est_total,1);
935  PCU_Add_Doubles(&u_norm_total,1);
936 
937  total_error = sqrt(err_est_total);
938  u_norm_total = sqrt(u_norm_total);
939  rel_err_total = total_error/u_norm_total;
940 
941  if(comm_rank==0){
942  std::cerr<<std::setprecision(10)<<std::endl;
943  std::cerr<<"Error estimate "<<total_error<<std::endl;
944  std::cerr<<"Error density maximum "<<errRho_max<<std::endl;
945  std::cerr<<"U_norm_total "<<u_norm_total<<std::endl;
946  }
947 
948  if(logging_config=="errorOnly"){ //feature to just look at the error fields without adapting the mesh
949  if(comm_rank==0)
950  std::cout<<"outputting error field\n";
951  char namebuffer[20];
952  sprintf(namebuffer,"err_reg_%i",nEstimate);
953  apf::writeVtkFiles(namebuffer, m);
954  //target_error = total_error*2; //this is a hack to prevent adapting
955  THRESHOLD = total_error*2;
956  nEstimate++;
957  removeBCData();
958  }
959 
960  m->end(iter);
961  apf::destroyField(visc);
962  apf::destroyField(estimate);
963 
964  if(comm_rank==0)
965  std::cerr<<"It cleared the ERM function.\n";
966 }
967 
MeshAdaptPUMIDrvr::total_error
double total_error
Definition: MeshAdaptPUMI.h:141
MeshAdaptPUMI.h
MeshAdaptPUMIDrvr::get_local_error
void get_local_error(double &total_error)
This function aims to compute error at each element via an Element Residual Method.
Definition: ErrorResidualMethod.cpp:664
MeshAdaptPUMIDrvr::nsd
int nsd
Definition: MeshAdaptPUMI.h:96
s
Double s
Definition: Headers.h:84
getLHS
void getLHS(Mat &K, apf::NewArray< apf::DynamicVector > &shdrv, int nsd, double weight, double visc_val, int nshl)
Definition: ErrorResidualMethod.cpp:197
MeshAdaptPUMIDrvr::localNumber
int localNumber(apf::MeshEntity *e)
Definition: MeshConverter.cpp:105
MeshAdaptPUMIDrvr::getMPvalue
double getMPvalue(double field_val, double val_0, double val_1)
Function primarily used to get the VOF-weighted average of physical properties at a given point.
Definition: ErrorResidualMethod.cpp:35
getDotProduct
double getDotProduct(apf::Vector3 a, apf::Vector3 b)
Definition: ErrorResidualMethod.cpp:73
setErrorField
void setErrorField(apf::Field *estimate, Vec coef, apf::MeshEntity *ent, int nsd, int nshl)
Definition: ErrorResidualMethod.cpp:601
MeshAdaptPUMIDrvr::getViscosityField
apf::Field * getViscosityField(apf::Field *voff)
Function used to derive a viscosity field from a VOF field.
Definition: ErrorResidualMethod.cpp:579
MeshAdaptPUMIDrvr::rel_err_total
double rel_err_total
Definition: MeshAdaptPUMI.h:143
MeshAdaptPUMIDrvr::nEstimate
int nEstimate
Definition: MeshAdaptPUMI.h:95
isInTet
bool isInTet(apf::Mesh *mesh, apf::MeshEntity *ent, apf::Vector3 pt)
Definition: ErrorResidualMethod.cpp:92
nu_0
double nu_0
Definition: ErrorResidualMethod.cpp:22
nu_1
double nu_1
Definition: ErrorResidualMethod.cpp:22
MeshAdaptPUMIDrvr::integration_order
int integration_order
Definition: MeshAdaptPUMI.h:138
MeshAdaptPUMIDrvr::approximation_order
int approximation_order
Definition: MeshAdaptPUMI.h:137
getRHS
void getRHS(Vec &F, apf::NewArray< double > &shpval, apf::NewArray< apf::DynamicVector > &shdrv, apf::Vector3 vel_vect, apf::Matrix3x3 grad_vel, int nsd, double weight, int nshl, double visc_val, double density, apf::Vector3 grad_density, double pressure, double g[3])
Definition: ErrorResidualMethod.cpp:244
MeshAdaptPUMIDrvr::computeDiffusiveFlux
void computeDiffusiveFlux(apf::Mesh *m, apf::Field *voff, apf::Field *visc, apf::Field *pref, apf::Field *velf)
Function used to compute the diffusive flux at interelement boundaries and stores the values as tags ...
Definition: ErrorResidualMethod.cpp:297
c
Double c
Definition: Headers.h:54
MeshAdaptPUMIDrvr::errRho_max
double errRho_max
Definition: MeshAdaptPUMI.h:142
int_order
int int_order
Definition: ErrorResidualMethod.cpp:21
rho_1
double rho_1
Definition: ErrorResidualMethod.cpp:22
MeshAdaptPUMIDrvr::BCtag
apf::MeshTag * BCtag
Definition: MeshAdaptPUMI.h:130
getProps
void getProps(double *rho, double *nu)
Definition: ErrorResidualMethod.cpp:25
a_kl
double a_kl
Definition: ErrorResidualMethod.cpp:23
MeshAdaptPUMIDrvr::fluxtag
apf::MeshTag * fluxtag[4]
Definition: MeshAdaptPUMI.h:132
MeshAdaptPUMIDrvr::logging_config
std::string logging_config
Definition: MeshAdaptPUMI.h:116
approx_order
int approx_order
Definition: ErrorResidualMethod.cpp:20
MeshAdaptPUMIDrvr::getBoundaryFlux
void getBoundaryFlux(apf::Mesh *m, apf::MeshEntity *ent, double *endflux)
This function reads in the stored tags and computes the boundary flux according to the RHS formulatio...
Definition: ErrorResidualMethod.cpp:494
isInSimplex
bool isInSimplex(apf::Mesh *mesh, apf::MeshEntity *ent, apf::Vector3 pt, int dim)
Definition: ErrorResidualMethod.cpp:127
rho_0
double rho_0
Definition: ErrorResidualMethod.cpp:22
getFaceNormal
apf::Vector3 getFaceNormal(apf::Mesh *mesh, apf::MeshEntity *face)
Definition: ErrorResidualMethod.cpp:50
MeshAdaptPUMIDrvr::removeBCData
void removeBCData()
Function used to remove the BC tags that were created during the computeDiffusiveFlux() function.
Definition: ErrorResidualMethod.cpp:626