proteus  1.8.1
C/C++/Fortran libraries
SizeField.cpp
Go to the documentation of this file.
1 #include "MeshAdaptPUMI.h"
2 #include <apf.h>
3 #include <apfVector.h>
4 #include <apfMesh.h>
5 #include <apfShape.h>
6 #include <apfDynamicVector.h>
7 #include <apfCavityOp.h>
8 #include <string>
9 #include <iostream>
10 #include <sstream>
11 #include <PCU.h>
12 #include <samElementCount.h>
13 #include <queue>
14 #include <algorithm>
15 
16 static void SmoothField(apf::Field *f);
17 void gradeAnisoMesh(apf::Mesh* m,double gradingFactor);
18 void gradeAspectRatio(apf::Mesh* m, int idx, double gradingFactor);
19 
20 /* Based on the distance from the interface epsilon can be controlled to determine
21  thickness of refinement near the interface */
22 static double isotropicFormula(double phi, double dphi, double verr, double hmin, double hmax, double phi_s = 0, double epsFact = 0)
23 {
24  double size;
25  double dphi_size_factor;
26  double v_size_factor;
27  if (fabs(phi_s) < (epsFact*7.5) * hmin)
28  return hmin;
29  else
30  return hmax;
31 }
32 
33 static void setSizeField(apf::Mesh2 *m,apf::MeshEntity *vertex,double h,apf::MeshTag *marker,apf::Field* sizeField)
34 //helper function for banded interface to facilitate with setting the proper mesh size and parallel communication
35 {
36  int isMarked=0;
37  if(m->hasTag(vertex,marker))
38  isMarked=1;
39  double h_new;
40  if(isMarked)
41  h_new = std::min(h,apf::getScalar(sizeField,vertex,0));
42  else
43  {
44  h_new = h;
45  int newMark = 1;
46  m->setIntTag(vertex,marker,&newMark);
47  }
48  apf::setScalar(sizeField,vertex,0,h_new);
49 
50  //Parallel Communication with owning copy
51  if(!m->isOwned(vertex))
52  {
53  apf::Copies remotes;
54  m->getRemotes(vertex,remotes);
55  int owningPart=m->getOwner(vertex);
56  PCU_COMM_PACK(owningPart, remotes[owningPart]);
57  PCU_COMM_PACK(owningPart, h_new);
58  }
59 }
60 
62 {
63  freeField(size_iso);
64  size_iso = apf::createLagrangeField(m, "proteus_size", apf::SCALAR, 1);
65 
66  apf::MeshIterator *it = m->begin(0);
67  apf::MeshEntity* ent;
68  while ((ent = m->iterate(it)))
69  {
70  int modelTag = m->getModelTag(m->toModel(ent));
71  //std::cout<<"This is the model tag "<<modelTag<<std::endl;
72  double sizeDesired;
73  if(modelTag==123)
74  sizeDesired=hmin;
75  else
76  sizeDesired=hmax;
77  apf::setScalar(size_iso,ent,0,sizeDesired);
78  }
79  m->end(it);
80  gradeMesh(1.5);
81 }
82 
83 
85 //Implementation of banded interface, edge intersection algorithm
86 //If mesh edge intersects the 0 level-set, then the adjacent edges need to be refined
87 {
88  if(m->findField("interfaceBand"))
89  apf::destroyField(m->findField("interfaceBand"));
90 
91  apf::Field* interfaceBand = apf::createLagrangeField(m, "interfaceBand", apf::SCALAR, 1);
92  apf::Field *phif = m->findField("phi");
93  assert(phif);
94 
95  apf::MeshTag* vertexMarker = m->createIntTag("vertexMarker",1);
96  apf::MeshIterator *it = m->begin(1);
97  apf::MeshEntity *edge;
98 
99  //double L_band = (numAdaptSteps+N_interface_band)*hPhi;
100 
101  PCU_Comm_Begin();
102  while ((edge = m->iterate(it)))
103  {
104  apf::Adjacent edge_adjVerts;
105  m->getAdjacent(edge,0,edge_adjVerts);
106  apf::MeshEntity *vertex1 = edge_adjVerts[0];
107  apf::MeshEntity *vertex2 = edge_adjVerts[1];
108  double phi1 = apf::getScalar(phif,vertex1,0);
109  double phi2 = apf::getScalar(phif,vertex2,0);
110  int caseNumber = 1;
111  if(std::fabs(phi1)>L_band)
112  caseNumber++;
113  if(std::fabs(phi2)>L_band)
114  caseNumber++;
115 
116  if(caseNumber==1 || caseNumber == 2)
117  {
118  setSizeField(m,vertex1,hPhi,vertexMarker,interfaceBand);
119  setSizeField(m,vertex2,hPhi,vertexMarker,interfaceBand);
120  }
121  else
122  {
123  if (phi1*phi2 <0)
124  {
125  setSizeField(m,vertex1,hPhi,vertexMarker,interfaceBand);
126  setSizeField(m,vertex2,hPhi,vertexMarker,interfaceBand);
127  }
128  else
129  {
130  setSizeField(m,vertex1,hmax,vertexMarker,interfaceBand);
131  setSizeField(m,vertex2,hmax,vertexMarker,interfaceBand);
132  }
133  }
134 
135  }//end while
136 
137  PCU_Comm_Send();
138 
139  //Take minimum between received value and current value
140  apf::MeshEntity *ent;
141  double h_received;
142  while(PCU_Comm_Receive())
143  {
144  //Note: the only receiving entities should be owning copies
145  PCU_COMM_UNPACK(ent);
146  PCU_COMM_UNPACK(h_received);
147  //take minimum of received values
148  double h_current = apf::getScalar(interfaceBand,ent,0);
149  double h_final = std::min(h_current,h_received);
150  apf::setScalar(interfaceBand,ent,0,h_final);
151  }
152 
153  //Synchronization has all remote copies track the owning copy value
154  apf::synchronize(interfaceBand);
155  m->end(it);
156 
157  m->destroyTag(vertexMarker);
158 
159  //add to queue
160  sizeFieldList.push(interfaceBand);
161  return 0;
162 }
163 
164 int intersectsInterface(apf::MeshEntity* edge, apf::Field* levelSet)
165 {
166  apf::Mesh* m = apf::getMesh(levelSet);
167  apf::Adjacent edge_adjVerts;
168  m->getAdjacent(edge,0,edge_adjVerts);
169  apf::MeshEntity *vertex1 = edge_adjVerts[0];
170  apf::MeshEntity *vertex2 = edge_adjVerts[1];
171  double phi1 = apf::getScalar(levelSet,vertex1,0);
172  double phi2 = apf::getScalar(levelSet,vertex2,0);
173  int doesIntersect = 0;
174  if(phi1*phi2 < 0) //implies different signs and therefore intersects interface
175  doesIntersect = 1;
176  return doesIntersect;
177 }
178 
179 //Struct definition
181  apf::MeshEntity* vertex;
182  apf::Vector3 actualPosition;
183  double direction;
184  double L_local;
185  apf::MeshTag* trackerTag;
186  //const char* tagName;
187  int edgeID;
189 };
190 
191 int checkForPropagation(apf::Mesh* m, edgeWalkerInfo inputObject)
192 {
193  apf::MeshEntity* vert = inputObject.vertex;
194  apf::Vector3 actualPosition = inputObject.actualPosition;
195  double L_local = inputObject.L_local;
196  double direction = inputObject.direction;
197 
198  apf::MeshTag* vertexMaxTraverse = m->findTag("maximumTraversal");
199  apf::Field* predictInterfaceBand = m->findField("predictInterfaceBand");
200  apf::Field* levelSet = m->findField("phi");
201 
202  apf::Vector3 pt_vert;
203  m->getPoint(vert,0,pt_vert);
204  apf::Vector3 difference_vect = pt_vert-actualPosition;
205 
206  //check if vertex needs to be added to queue based on traversal distance
207  int dontContinue = 0;
208  if(m->hasTag(vert,vertexMaxTraverse))
209  {
210  double traversalDistance;
211  m->getDoubleTag(vert,vertexMaxTraverse,&traversalDistance);
212  if((L_local-difference_vect.getLength()) < traversalDistance*1.01) //has to be 1% higher
213  dontContinue=1;
214  }
215 
216  //directionality
217  double phiCurrent = apf::getScalar(levelSet,vert,0);
218  if((difference_vect.getLength() > L_local) || dontContinue || (phiCurrent*direction<=0))
219  return 0;
220  else
221  return 1;
222 }
223 
224 int BFS_propagation(apf::Mesh* m, std::queue<edgeWalkerInfo> &markedVertices)
225 {
226 //objects in the queue are assumed to be checked and needing to be modified
227 //the adjacencies are checked before adding into the queue
228 //for parallelism, what's important is to stop a constant back-and-forth communication on shared vertices; we do this by considering communication when we add an entry into the queue
229 
230  //get the latest object
231  edgeWalkerInfo inputObject = markedVertices.front();
232  markedVertices.pop();
233 
234  //set the variables from the inputObject
235  apf::MeshEntity* vert = inputObject.vertex;
236  apf::Vector3 actualPosition = inputObject.actualPosition;
237  double L_local = inputObject.L_local;
238  double direction = inputObject.direction;
239 
240  //get necessary fields
241  apf::MeshTag* vertexMaxTraverse = m->findTag("maximumTraversal");
242  apf::Field* predictInterfaceBand = m->findField("predictInterfaceBand");
243  apf::Field* levelSet = m->findField("phi");
244 
245  int needsParallel=0;
246 
247  apf::Adjacent vertex_adjVerts;
248  apf::getBridgeAdjacent(m,vert,1,0,vertex_adjVerts);
249  for(int i=0;i<vertex_adjVerts.getSize();i++)
250  {
251 
252  apf::MeshEntity* newVert = vertex_adjVerts[i];
253  inputObject.vertex=newVert;
254 
255  if(checkForPropagation(m,inputObject))
256  {
257 
258  //set new traversal distance
259  apf::Vector3 pt_vert;
260  m->getPoint(newVert,0,pt_vert);
261  apf::Vector3 difference_vect = pt_vert-actualPosition;
262 
263  double traversalDistance = L_local-difference_vect.getLength();
264  m->setDoubleTag(newVert,vertexMaxTraverse,&traversalDistance);
265 
266  //set size
267  apf::setScalar(predictInterfaceBand,newVert,0,apf::getScalar(predictInterfaceBand,vert,0));
268  //add vertex to queue
269  inputObject.vertex = newVert;
270  markedVertices.push(inputObject);
271 
272  if(m->isShared(newVert))
273  {
274  int initialRank = PCU_Comm_Self();
275  double desiredSize = apf::getScalar(predictInterfaceBand,vert,0);
276  apf::Copies remotes;
277  m->getRemotes(newVert,remotes);
278  for(apf::Copies::iterator iter=remotes.begin(); iter!=remotes.end();++iter)
279  {
280  PCU_COMM_PACK(iter->first, iter->second);
281  PCU_COMM_PACK(iter->first, L_local);
282  PCU_COMM_PACK(iter->first, actualPosition);
283  PCU_COMM_PACK(iter->first, direction);
284  PCU_COMM_PACK(iter->first, initialRank);
285  PCU_COMM_PACK(iter->first, inputObject.edgeID);
286  PCU_COMM_PACK(iter->first, desiredSize);
287  }
288  needsParallel++;
289  }
290  }
291  } //end of adjacent
292 
293  return needsParallel;
294 }
295 
296 
298 //compute Lband
299 //edge walk
300 {
301  apf::Field* interfaceBand = m->findField("interfaceBand");
302  apf::Field* velocity = m->findField("velocity");
303  apf::Field* levelSet = m->findField("phi");
304 
305  //get gradient field
306  apf::Field *gradphi = apf::recoverGradientByVolume(levelSet);
307 
308  apf::Field* predictInterfaceBand = apf::createLagrangeField(m,"predictInterfaceBand",apf::SCALAR,1);
309  apf::copyData(predictInterfaceBand,interfaceBand);
310 
311  //edge-walk to predict
312  apf::MeshTag* vertexMaxTraverse = m->createDoubleTag("maximumTraversal",1); //define tag field for each vertex to store maximum distance that will be travelled
313 
314  apf::MeshEntity* edge;
315  apf::MeshIterator* it = m->begin(1);
316 
317  std::queue <edgeWalkerInfo> markedVertices;
318 
319  PCU_Comm_Begin();
320  while( (edge = m->iterate(it)) )
321  {
322  if( intersectsInterface(edge,levelSet))
323  {
324  //get the parameterized position of interface along edge
325  apf::Adjacent edge_adjVerts;
326  m->getAdjacent(edge,0,edge_adjVerts);
327  apf::MeshEntity *vertex1 = edge_adjVerts[0];
328  apf::MeshEntity *vertex2 = edge_adjVerts[1];
329  double phi1 = apf::getScalar(levelSet,vertex1,0);
330  double phi2 = apf::getScalar(levelSet,vertex2,0);
331  apf::Vector3 pt_1, pt_2;
332  m->getPoint(vertex1,0,pt_1);
333  m->getPoint(vertex2,0,pt_2);
334  apf::Vector3 edgeVector = pt_2-pt_1;
335  double edgeLength = apf::measure(m,edge);
336  double zeroPosition = 2*(-phi1/(phi2-phi1))-1.0; //parametric position of interface on the edge
337  double relativePosition = -phi1/(phi2-phi1); //same as zeroPosition but in interval of [0,1]
338  apf::Vector3 actualPosition = (pt_2-pt_1)*relativePosition + pt_1;
339 
340  apf::Vector3 edgePoint(zeroPosition,0.0,0.0);
341  apf::Element* phiElem = apf::createElement(levelSet,edge);
342  apf::Element* gradPhiElem = apf::createElement(gradphi,edge);
343  apf::Element* velocityElem = apf::createElement(velocity,edge);
344  apf::Vector3 localVelocity;
345  apf::getVector(velocityElem,edgePoint,localVelocity);
346  apf::Vector3 localInterfaceNormal;
347  apf::getVector(gradPhiElem,edgePoint,localInterfaceNormal);
348  apf::destroyElement(phiElem);
349  apf::destroyElement(velocityElem);
350  apf::destroyElement(gradPhiElem);
351 
352  //get L_local
353  double L_local = localVelocity.getLength()*numAdaptSteps*delta_T;
354 
355  L_local += (N_interface_band)*hPhi; //add blending region
356 
357  //get direction, multiply this with levelSet value to determine if in same direction
358  double signValue = localVelocity*localInterfaceNormal;
359 
360  //find adjacent vertices and their adjacent edges
361  for(int i=0; i<edge_adjVerts.getSize();i++)
362  {
363  edgeWalkerInfo inputObject;
364  inputObject.vertex = edge_adjVerts[i];
365  inputObject.actualPosition = actualPosition;
366  inputObject.direction = signValue;
367  inputObject.L_local = L_local;
368  inputObject.edgeID = localNumber(edge);
369  inputObject.initialRank = PCU_Comm_Self();
370  if(checkForPropagation(m,inputObject))
371  {
372  markedVertices.push(inputObject);
373 
374  //check for parallel
375  if(m->isShared(inputObject.vertex))
376  {
377  int initialRank = PCU_Comm_Self();
378  double desiredSize = apf::getScalar(predictInterfaceBand,inputObject.vertex,0);
379  apf::Copies remotes;
380  m->getRemotes(inputObject.vertex,remotes);
381  for(apf::Copies::iterator iter=remotes.begin(); iter!=remotes.end();++iter)
382  {
383  PCU_COMM_PACK(iter->first, iter->second);
384  PCU_COMM_PACK(iter->first, L_local);
385  PCU_COMM_PACK(iter->first, actualPosition);
386  PCU_COMM_PACK(iter->first, inputObject.direction);
387  PCU_COMM_PACK(iter->first, initialRank);
388  PCU_COMM_PACK(iter->first, inputObject.edgeID);
389  PCU_COMM_PACK(iter->first, desiredSize);
390  }
391  }
392 
393  }
394  }
395  } //end if interface edge
396  }
397  m->end(it);
398 
399  //The following parallel code is just for the initialization step
400  //There might be a way to put all of this under a function as this is repeated code later on
401  PCU_Comm_Send();
402  while(PCU_Comm_Receive())
403  {
404  apf::MeshEntity* vertex;
405  double L_local;
406  apf::Vector3 actualPosition;
407  double direction;
408  int initialRank;
409  int edgeID;
410  double desiredSize;
411  PCU_COMM_UNPACK(vertex);
412  PCU_COMM_UNPACK(L_local);
413  PCU_COMM_UNPACK(actualPosition);
414  PCU_COMM_UNPACK(direction);
415  PCU_COMM_UNPACK(initialRank);
416  PCU_COMM_UNPACK(edgeID);
417  PCU_COMM_UNPACK(desiredSize);
418 
419  edgeWalkerInfo inputObject;
420  inputObject.vertex = vertex;
421  inputObject.L_local = L_local;
422  inputObject.actualPosition = actualPosition;
423  inputObject.direction = direction;
424  inputObject.edgeID = edgeID;
425  inputObject.initialRank=initialRank;
426 
427  //ensures that the size is the same across parts
428  if(desiredSize < apf::getScalar(predictInterfaceBand,vertex,0))
429  apf::setScalar(predictInterfaceBand,vertex,0,desiredSize);
430 
431  if(checkForPropagation(m,inputObject))
432  {
433  apf::Vector3 pt_vert;
434  m->getPoint(vertex,0,pt_vert);
435  apf::Vector3 difference_vect = pt_vert-actualPosition;
436 
437  double traversalDistance = L_local-difference_vect.getLength();
438  m->setDoubleTag(vertex,vertexMaxTraverse,&traversalDistance);
439 
440  markedVertices.push(inputObject);
441  }
442  }
443 
444  //Parallel preparations
445  int needsParallel=1;
446 
447  while(needsParallel>0)
448  {
449  needsParallel=0;
450 
451  PCU_Comm_Begin();
452 
453  //Handle the queue
454  while(!markedVertices.empty())
455  {
456  needsParallel+=BFS_propagation(m,markedVertices);
457  }
458  PCU_Add_Ints(&needsParallel,1);
459 
460  PCU_Comm_Send();
461  while( PCU_Comm_Receive() )
462  {
463  apf::MeshEntity* vertex;
464  double L_local;
465  apf::Vector3 actualPosition;
466  double direction;
467  int initialRank;
468  int edgeID;
469  double desiredSize;
470  PCU_COMM_UNPACK(vertex);
471  PCU_COMM_UNPACK(L_local);
472  PCU_COMM_UNPACK(actualPosition);
473  PCU_COMM_UNPACK(direction);
474  PCU_COMM_UNPACK(initialRank);
475  PCU_COMM_UNPACK(edgeID);
476  PCU_COMM_UNPACK(desiredSize);
477 
478 
479  edgeWalkerInfo inputObject;
480  inputObject.vertex = vertex;
481  inputObject.L_local = L_local;
482  inputObject.actualPosition = actualPosition;
483  inputObject.direction = direction;
484  inputObject.edgeID = edgeID;
485  inputObject.initialRank=initialRank;
486 
487  //ensures that the size is the same across parts
488  if(desiredSize < apf::getScalar(predictInterfaceBand,vertex,0))
489  apf::setScalar(predictInterfaceBand,vertex,0,desiredSize);
490 
491  if(checkForPropagation(m,inputObject))
492  {
493  apf::Vector3 pt_vert;
494  m->getPoint(vertex,0,pt_vert);
495  apf::Vector3 difference_vect = pt_vert-actualPosition;
496 
497  double traversalDistance = L_local-difference_vect.getLength();
498  m->setDoubleTag(vertex,vertexMaxTraverse,&traversalDistance);
499 
500  markedVertices.push(inputObject);
501  }
502  }
503  }
504 
505  apf::copyData(interfaceBand,predictInterfaceBand);
506  apf::destroyField(predictInterfaceBand);
507  apf::destroyField(gradphi);
508 
509  //get all tags for destruction
510  apf::MeshIterator* tagIt = m->begin(0);
511  apf::MeshEntity* taggedVertex;
512  while( taggedVertex = m->iterate(tagIt) )
513  {
514  if(m->hasTag(taggedVertex,vertexMaxTraverse))
515  m->removeTag(taggedVertex,vertexMaxTraverse);
516  }
517  m->end(tagIt);
518  m->destroyTag(vertexMaxTraverse);
519 
520 }
521 
522 void MeshAdaptPUMIDrvr::isotropicIntersect()
523 {
524  freeField(size_iso);
525  if(m->findField("proteus_size"))
526  apf::destroyField(m->findField("proteus_size"));
527 
528  size_iso = apf::createFieldOn(m, "proteus_size", apf::SCALAR);
529 
530  apf::MeshEntity *vert;
531  apf::MeshIterator *it = m->begin(0);
532 
533  apf::Field *field = sizeFieldList.front();
534  apf::copyData(size_iso,field);
535  sizeFieldList.pop();
536  //apf::destroyField(field);
537  while(!sizeFieldList.empty())
538  {
539  field = sizeFieldList.front();
540  while(vert = m->iterate(it))
541  {
542  double value1 = apf::getScalar(size_iso,vert,0);
543  double value2 = apf::getScalar(field,vert,0);
544  double minValue = std::min(value1,value2);
545  apf::setScalar(size_iso,vert,0,minValue);
546  }
547  sizeFieldList.pop();
548  //apf::destroyField(field);
549  }
550  if(nAdapt < 2) //if first few adapts, need to make sure that gradation is low to allow interface to develop properly
551  gradeMesh(1.1);
552  else
554 }
555 
556 //taken from Dan's superconvergent patch recovery code
557 void MeshAdaptPUMIDrvr::averageToEntity(apf::Field *ef, apf::Field *vf,
558  apf::MeshEntity *ent)
559 //Function used to convert a region/element field into a vertex field via averaging.
560 //For each vertex, the average value of the adjacent elements is taken
561 //Input:
562 // ef is the element field
563 // ent is the target vertex
564 //Output:
565 // vf is the vertex field
566 {
567  apf::Mesh *m = apf::getMesh(ef);
568  apf::Adjacent elements;
569  m->getAdjacent(ent, m->getDimension(), elements);
570  double s = 0;
571  for (std::size_t i = 0; i < elements.getSize(); ++i)
572  s += apf::getScalar(ef, elements[i], 0);
573  s /= elements.getSize();
574  apf::setScalar(vf, ent, 0, s);
575  return;
576 }
577 
578 void MeshAdaptPUMIDrvr::minToEntity(apf::Field* ef, apf::Field* vf,
579  apf::MeshEntity* ent)
580 {
581  apf::Mesh* m = apf::getMesh(ef);
582  apf::Adjacent elements;
583  m->getAdjacent(ent, m->getDimension(), elements);
584  double s=0;
585  for (std::size_t i=0; i < elements.getSize(); ++i){
586  if(i==0)
587  s = apf::getScalar(ef, elements[i], 0);
588  else if(apf::getScalar(ef,elements[i],0) < s)
589  s= apf::getScalar(ef,elements[i],0);
590  }
591  apf::setScalar(vf, ent, 0, s);
592  return;
593 }
594 
595 
596 void MeshAdaptPUMIDrvr::volumeAverageToEntity(apf::Field *ef, apf::Field *vf,
597  apf::MeshEntity *ent)
598 //Serves the same purpose as averageToEntity but considers a volume-weighted average
599 {
600  apf::Mesh *m = apf::getMesh(ef);
601  apf::Adjacent elements;
602  apf::MeshElement *testElement;
603  m->getAdjacent(ent, m->getDimension(), elements);
604  double s = 0;
605  double VolumeTotal = 0;
606  for (std::size_t i = 0; i < elements.getSize(); ++i)
607  {
608  testElement = apf::createMeshElement(m, elements[i]);
609  s += apf::getScalar(ef, elements[i], 0)*apf::measure(testElement);
610  VolumeTotal += apf::measure(testElement);
611  apf::destroyMeshElement(testElement);
612  }
613  s /= VolumeTotal;
614  apf::setScalar(vf, ent, 0, s);
615  return;
616 }
617 
618 void errorAverageToEntity(apf::Field *ef, apf::Field *vf, apf::Field* err, apf::MeshEntity *ent)
619 //Serves the same purpose as averageToEntity but considers a error-weighted average
620 {
621  apf::Mesh *m = apf::getMesh(ef);
622  apf::Adjacent elements;
623  m->getAdjacent(ent, m->getDimension(), elements);
624  double s = 0;
625  double errorTotal = 0;
626  for (std::size_t i = 0; i < elements.getSize(); ++i)
627  {
628  s += apf::getScalar(ef, elements[i], 0)*apf::getScalar(err,elements[i],0);
629  errorTotal += apf::getScalar(err,elements[i],0);
630  }
631  s /= errorTotal;
632 /*
633  if (comm_rank == 0)
634  {
635  std::cout << "What is s final? " << s << std::endl;
636  }
637 */
638  apf::setScalar(vf, ent, 0, s);
639  return;
640 }
641 
642 
643 static apf::Field *extractSpeed(apf::Field *velocity)
644 //Function used to convert the velocity field into a speed field
645 {
646  apf::Mesh *m = apf::getMesh(velocity);
647  apf::Field *speedF = apf::createLagrangeField(m, "proteus_speed", apf::SCALAR, 1);
648  apf::MeshIterator *it = m->begin(0);
649  apf::MeshEntity *v;
650  apf::Vector3 vel_vect;
651  while ((v = m->iterate(it)))
652  {
653  apf::getVector(velocity, v, 0, vel_vect);
654  double speed = vel_vect.getLength();
655  apf::setScalar(speedF, v, 0, speed);
656  }
657  m->end(it);
658  return speedF;
659 }
660 
661 static apf::Matrix3x3 hessianFormula(apf::Matrix3x3 const &g2phi)
662 //Function used to output a symmetric hessian
663 {
664  apf::Matrix3x3 g2phit = apf::transpose(g2phi);
665  return (g2phi + g2phit) / 2;
666 }
667 
668 static apf::Field *computeHessianField(apf::Field *grad2phi)
669 //Function that writes a hessian field
670 {
671  apf::Mesh *m = apf::getMesh(grad2phi);
672  apf::Field *hessf = createLagrangeField(m, "proteus_hess", apf::MATRIX, 1);
673  apf::MeshIterator *it = m->begin(0);
674  apf::MeshEntity *v;
675  while ((v = m->iterate(it)))
676  {
677  apf::Matrix3x3 g2phi;
678  apf::getMatrix(grad2phi, v, 0, g2phi);
679  apf::Matrix3x3 hess = hessianFormula(g2phi);
680  apf::setMatrix(hessf, v, 0, hess);
681  }
682  m->end(it);
683  return hessf;
684 }
685 
686 static apf::Field *computeMetricField(apf::Field *gradphi, apf::Field *grad2phi, apf::Field *size_iso, double eps_u)
687 //Function used to generate a field of metric tensors at vertices. It is meant to be an umbrella function that can compute Hessians too.
688 //Currently needs development.
689 {
690  apf::Mesh *m = apf::getMesh(grad2phi);
691  apf::Field *metricf = createLagrangeField(m, "proteus_metric", apf::MATRIX, 1);
692  apf::MeshIterator *it = m->begin(0);
693  apf::MeshEntity *v;
694  while ((v = m->iterate(it)))
695  {
696  apf::Matrix3x3 g2phi;
697  apf::getMatrix(grad2phi, v, 0, g2phi);
698 /*
699  apf::Vector3 gphi;
700  apf::getVector(gradphi, v, 0, gphi);
701  apf::Matrix3x3 gphigphit(gphi[0] * gphi[0], gphi[0] * gphi[1], gphi[0] * gphi[2],
702  gphi[0] * gphi[1], gphi[1] * gphi[1], gphi[1] * gphi[2],
703  gphi[0] * gphi[2], gphi[1] * gphi[2], gphi[2] * gphi[2]);
704 */
705  apf::Matrix3x3 hess = hessianFormula(g2phi);
706  apf::Matrix3x3 metric = hess;
707  //apf::Matrix3x3 metric = gphigphit/(apf::getScalar(size_iso,v,0)*apf::getScalar(size_iso,v,0))+ hess/eps_u;
708  //apf::Matrix3x3 metric(1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0);
709  apf::setMatrix(metricf, v, 0, metric);
710  }
711  m->end(it);
712  return metricf;
713 }
714 
715 // Gaussian, Mean and principal curvatures based on Hessian and gradient of phi.
716 static void curveFormula(apf::Matrix3x3 const &h, apf::Vector3 const &g,
717  apf::Vector3 &curve)
718 {
719  double a = (h[1][1] + h[2][2]) * g[0] * g[0] + (h[0][0] + h[2][2]) * g[1] * g[1] + (h[0][0] + h[1][1]) * g[2] * g[2];
720 
721  double b = g[0] * g[1] * h[0][1] + g[0] * g[2] * h[0][2] + g[1] * g[2] * h[1][2];
722 
723  double Km = (a - 2 * b) / pow(g * g, 1.5);
724 
725  double c = g[0] * g[0] * (h[1][1] * h[2][2] - h[1][2] * h[1][2]) + g[1] * g[1] * (h[0][0] * h[2][2] - h[0][2] * h[0][2]) + g[2] * g[2] * (h[0][0] * h[1][1] - h[0][1] * h[0][1]);
726 
727  double d = g[0] * g[1] * (h[0][2] * h[1][2] - h[0][1] * h[2][2]) + g[1] * g[2] * (h[0][1] * h[0][2] - h[1][2] * h[0][0]) + g[0] * g[2] * (h[0][1] * h[1][2] - h[0][2] * h[1][1]);
728 
729  double Kg = (c + 2 * d) / pow(g * g, 2);
730 
731  curve[0] = Km; //Mean curvature= (k_1+k_2)/2, Not used in normal direction
732  curve[1] = Km + sqrt(Km * Km - Kg); //k_1, First principal curvature (maximum curvature)
733  curve[2] = Km - sqrt(Km * Km - Kg); //k_2, Second principal curvature (minimum curvature)
734 }
735 
736 static apf::Field *getCurves(apf::Field *hessians, apf::Field *gradphi)
737 {
738  apf::Mesh *m = apf::getMesh(hessians);
739  apf::Field *curves;
740  curves = apf::createLagrangeField(m, "proteus_curves", apf::VECTOR, 1);
741  apf::MeshIterator *it = m->begin(0);
742  apf::MeshEntity *v;
743  while ((v = m->iterate(it)))
744  {
745  apf::Matrix3x3 hessian;
746  apf::getMatrix(hessians, v, 0, hessian);
747  apf::Vector3 gphi;
748  apf::getVector(gradphi, v, 0, gphi);
749  apf::Vector3 curve;
750  curveFormula(hessian, gphi, curve);
751  apf::setVector(curves, v, 0, curve);
752  }
753  m->end(it);
754  return curves;
755 }
756 
757 static void clamp(double &v, double min, double max)
758 //Function used to restrict the mesh size to user specified minimums and maximums
759 {
760  v = std::min(max, v);
761  v = std::max(min, v);
762 }
763 
764 static void clampField(apf::Field *field, double min, double max)
765 //Function that loops through a field and clamps the values
766 {
767  apf::Mesh *m = apf::getMesh(field);
768  int numcomps = apf::countComponents(field);
769  double components[numcomps];
770  apf::MeshEntity *v;
771  apf::MeshIterator *it = m->begin(0);
772  while ((v = m->iterate(it)))
773  {
774  //double tempValue = apf::getScalar(field,v,0);
775  apf::getComponents(field, v, 0, &components[0]);
776  for (int i = 0; i < numcomps; i++)
777  clamp(components[i], min, max);
778  //apf::setScalar(field,v,0,tempValue);
779  apf::setComponents(field, v, 0, &components[0]);
780  }
781  m->end(it);
782 }
783 
784 static void scaleFormula(double phi, double hmin, double hmax,
785  int adapt_step,
786  apf::Vector3 const &curves,
787  apf::Vector3 &scale)
788 //Function used to set the size scale vector for the interface size field configuration
789 {
790  double epsilon = 7.0 * hmin;
791  if (fabs(phi) < epsilon)
792  {
793  scale[0] = hmin;
794  scale[1] = sqrt(0.002 / fabs(curves[1]));
795  scale[2] = sqrt(0.002 / fabs(curves[2]));
796  }
797  else if (fabs(phi) < 4 * epsilon)
798  {
799  scale[0] = 2 * hmin;
800  scale[1] = 2 * sqrt(0.002 / fabs(curves[1]));
801  scale[2] = 2 * sqrt(0.002 / fabs(curves[2]));
802  }
803  else
804  {
805  scale = apf::Vector3(1, 1, 1) * hmax;
806  }
807 
808  //for (int i = 0; i < 3; ++i)
809  // clamp(scale[i], hmin, hmax);
810 }
811 
812 static void scaleFormulaERM(double phi, double hmin, double hmax, double h_dest,
813  apf::Vector3 const &curves,
814  double lambda[3], double eps_u, apf::Vector3 &scale,int nsd,double maxAspect)
815 //Function used to set the size scale vector for the anisotropic ERM size field configuration
816 //Inputs:
817 // phi is is the distance to the interface
818 // hmin is the minimum mesh size
819 // hmax is the maximum mesh size
820 // h_dest is the computed new mesh size
821 // curves is a vector that denotes the curvature of the interface
822 // lambda is the ordered set of eigenvalues from an eigendecomposition of the metric tensor
823 // eps_u is a tolerance for the distance away from the interface
824 //Output:
825 // scale is the mesh size in each direction for a vertex
826 {
827 /*
828  double epsilon = 7.0 * hmin;
829  double lambdamin = 1.0 / (hmin * hmin);
830  if (lambda[1] < 1e-10)
831  {
832  lambda[1] = lambdamin;
833  lambda[2] = lambdamin;
834  }
835  if (lambda[2] < 1e-10)
836  {
837  lambda[2] = lambdamin;
838  }
839 */
841 
842 /*
843  scale[0] = h_dest*pow(lambda[1]/lambda[0],0.25);
844  scale[1] = sqrt(lambda[0]/lambda[1])*scale[0];
845  scale[2] = 1.0;
846 */
847 /*
848  scale[0] = h_dest*pow(lambda[1]/lambda[0],0.25)*pow(3,0.25)*0.5;
849  scale[1] = sqrt(lambda[0]/lambda[1])*scale[0];
850  scale[2] = 1.0;
851 */
852 /*
853  if(nsd == 2){
854  scale[0] = h_dest;
855  scale[1] = sqrt(lambda[0]/lambda[1])*scale[0];
856  scale[2] = 1.0;
857  }
858  else{
859  scale[0] = h_dest;
860  scale[1] = sqrt(lambda[0]/lambda[1])*scale[0];
861  scale[2] = sqrt(lambda[0]/lambda[2])*scale[0];
862  }
863 */
864 
865 //3D
866 
867  if(nsd == 2){
868  scale[0] = h_dest * pow((lambda[1] ) / (lambda[0]), 1.0 / 4.0);
869  scale[1] = sqrt(lambda[0] / lambda[1]) * scale[0];
870  scale[2] = 1.0;
871  }
872  else{
873 /*
874  scale[0] = h_dest * pow((lambda[1] * lambda[2]) / (lambda[0] * lambda[0]), 1.0 / 6.0);
875  scale[1] = sqrt(lambda[0] / lambda[1]) * scale[0];
876  scale[2] = sqrt(lambda[0] / lambda[2]) * scale[0];
877 */
878  scale[0] = h_dest;
879  scale[1] = sqrt(lambda[0] / lambda[1]) * scale[0];
880  scale[2] = sqrt(lambda[0] / lambda[2]) * scale[0];
881  if(scale[1]/scale[0] > maxAspect)
882  scale[1] = maxAspect*scale[0];
883  if(scale[2]/scale[0] > maxAspect)
884  scale[2] = maxAspect*scale[0];
885  if(scale[1]/scale[0] > maxAspect || scale[2]/scale[0] > maxAspect){
886  logEvent("Scales reached maximum aspect ratio",4);
887  }
888  }
889  //*/
890  /*
891  if(fabs(phi)<epsilon){
892  scale[0] = h_dest;
893  scale[1] = sqrt(eps_u/fabs(curves[2]));
894  scale[2] = sqrt(eps_u/fabs(curves[1]));
895  }
896  else
897  scale = apf::Vector3(1,1,1) * h_dest;
898 */
899 }
900 
901 static apf::Field *getSizeScales(apf::Field *phif, apf::Field *curves,
902  double hmin, double hmax, int adapt_step)
903 //Function that gets size scales in the interface size field configuration
904 {
905  apf::Mesh *m = apf::getMesh(phif);
906  apf::Field *scales;
907  scales = apf::createLagrangeField(m, "proteus_size_scale", apf::VECTOR, 1);
908  apf::MeshIterator *it = m->begin(0);
909  apf::MeshEntity *v;
910  while ((v = m->iterate(it)))
911  {
912  double phi = apf::getScalar(phif, v, 0);
913  apf::Vector3 curve;
914  apf::getVector(curves, v, 0, curve);
915  apf::Vector3 scale;
916  scaleFormula(phi, hmin, hmax, adapt_step, curve, scale);
917  apf::setVector(scales, v, 0, scale);
918  }
919  m->end(it);
920  return scales;
921 }
922 
924 {
925  apf::Vector3 v;
926  double wm;
927  bool operator<(const SortingStruct &other) const
928  {
929  return wm < other.wm;
930  }
931 };
932 
933 static apf::Field *getSizeFrames(apf::Field *hessians, apf::Field *gradphi)
934 //Function that gets size frames, or the metric tensor, in the interface size field configuration
935 {
936  apf::Mesh *m = apf::getMesh(gradphi);
937  apf::Field *frames;
938  frames = apf::createLagrangeField(m, "proteus_size_frame", apf::MATRIX, 1);
939  apf::MeshIterator *it = m->begin(0);
940  apf::MeshEntity *v;
941  while ((v = m->iterate(it)))
942  {
943  apf::Vector3 gphi;
944  apf::getVector(gradphi, v, 0, gphi);
945  apf::Vector3 dir;
946  if (gphi.getLength() > 1e-16)
947  dir = gphi.normalize();
948  else
949  dir = apf::Vector3(1, 0, 0);
950  apf::Matrix3x3 hessian;
951  apf::getMatrix(hessians, v, 0, hessian);
952  apf::Vector3 eigenVectors[3];
953  double eigenValues[3];
954  apf::eigen(hessian, eigenVectors, eigenValues);
955  SortingStruct ssa[3];
956  for (int i = 0; i < 3; ++i)
957  {
958  ssa[i].v = eigenVectors[i];
959  ssa[i].wm = std::fabs(eigenValues[i]);
960  }
961  std::sort(ssa, ssa + 3);
962  assert(ssa[2].wm >= ssa[1].wm);
963  assert(ssa[1].wm >= ssa[0].wm);
964  double firstEigenvalue = ssa[2].wm;
965  apf::Matrix3x3 frame;
966  frame[0] = dir;
967  if (firstEigenvalue > 1e-16)
968  {
969  apf::Vector3 firstEigenvector = ssa[2].v;
970  frame[1] = apf::reject(firstEigenvector, dir);
971  frame[2] = apf::cross(frame[0], frame[1]);
972  if (frame[2].getLength() < 1e-16)
973  frame = apf::getFrame(dir);
974  }
975  else
976  frame = apf::getFrame(dir);
977  for (int i = 0; i < 3; ++i)
978  frame[i] = frame[i].normalize();
979  frame = apf::transpose(frame);
980  apf::setMatrix(frames, v, 0, frame);
981  }
982  m->end(it);
983  return frames;
984 }
985 
986 static apf::Field *getERMSizeFrames(apf::Field *hessians, apf::Field *gradphi, apf::Field *frame_comps[3])
987 //Function that gets size frames, or the metric tensor, in the interface size field configuration
988 {
989  apf::Mesh *m = apf::getMesh(gradphi);
990  apf::Field *frames;
991  frames = m->findField("proteus_size_frame");
992  apf::MeshIterator *it = m->begin(0);
993  apf::MeshEntity *v;
994  while ((v = m->iterate(it)))
995  {
996  apf::Matrix3x3 frame(1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0);
997  apf::Vector3 gphi;
998  apf::getVector(gradphi, v, 0, gphi);
999 /*
1000  apf::Vector3 dir;
1001  if (gphi.getLength() > 1e-16)
1002  dir = gphi.normalize();
1003  else
1004  dir = apf::Vector3(1, 0, 0);
1005 */
1006 
1007  //get eigen values and eigenvectors from hessian
1008  apf::Matrix3x3 hessian;
1009  apf::getMatrix(hessians, v, 0, hessian);
1010  apf::Vector3 eigenVectors[3];
1011  double eigenValues[3];
1012  apf::eigen(hessian, eigenVectors, eigenValues);
1013  SortingStruct ssa[3];
1014  for (int i = 0; i < 3; ++i)
1015  {
1016  ssa[i].v = eigenVectors[i];
1017  ssa[i].wm = std::fabs(eigenValues[i]);
1018  }
1019 
1020  //sort eigenvalues and eigenvectors
1021  std::sort(ssa, ssa + 3);
1022  assert(ssa[2].wm >= ssa[1].wm);
1023  assert(ssa[1].wm >= ssa[0].wm);
1024  double firstEigenvalue = ssa[2].wm;
1025  assert(firstEigenvalue > 1e-12);
1026  //frame[0] = dir;
1027  //if (firstEigenvalue > 1e-16)
1028  //{
1029  frame[0] = ssa[2].v;
1030  frame[1] = ssa[1].v;
1031  frame[2] = ssa[0].v;
1032  //}
1033  //else
1034  // frame = apf::getFrame(dir);
1035 
1036 /*
1037  apf::Vector3 test(1.0,0.0,0.0);
1038  frame = apf::getFrame(test);
1039  apf::setMatrix(frames,v,0,frame);
1040 */
1041 
1042  //normalize eigenvectors
1043  for (int i = 0; i < 3; ++i)
1044  frame[i] = frame[i].normalize();
1045  frame = apf::transpose(frame);
1046  apf::setMatrix(frames, v, 0, frame);
1047  apf::setVector(frame_comps[0], v, 0, frame[0]);
1048  apf::setVector(frame_comps[1], v, 0, frame[1]);
1049  apf::setVector(frame_comps[2], v, 0, frame[2]);
1050  }
1051  m->end(it);
1052  return frames;
1053 }
1054 
1056 //High level function that obtains the size scales and the size frames for anistropic interface-based adapt
1057 {
1058  apf::Field *phif = m->findField("phi");
1059  assert(phif);
1060  apf::Field *gradphi = apf::recoverGradientByVolume(phif);
1061  apf::Field *grad2phi = apf::recoverGradientByVolume(gradphi);
1062  apf::Field *hess = computeHessianField(grad2phi);
1063  apf::destroyField(grad2phi);
1064  apf::Field *curves = getCurves(hess, gradphi);
1065  freeField(size_scale);
1066 
1067  size_scale = getSizeScales(phif, curves, hmin, hmax, nAdapt);
1068  apf::destroyField(curves);
1069  freeField(size_frame);
1070  size_frame = getSizeFrames(hess, gradphi);
1071  apf::destroyField(hess);
1072 
1073  apf::destroyField(gradphi);
1074  for (int i = 0; i < 2; ++i)
1075  SmoothField(size_scale);
1076 
1077  return 0;
1078 }
1079 
1080 struct Smoother : public apf::CavityOp
1081 {
1082  Smoother(apf::Field *f) : apf::CavityOp(apf::getMesh(f))
1083  {
1084  field = f;
1085  int nc = apf::countComponents(f);
1086  newField = apf::createPackedField(mesh, "proteus_smooth_new", nc);
1087  sum.setSize(nc);
1088  value.setSize(nc);
1089  nApplied = 0;
1090  }
1092  {
1093  copyData(field, newField);
1094  apf::destroyField(newField);
1095  }
1096  virtual Outcome setEntity(apf::MeshEntity *e)
1097  {
1098  if (apf::hasEntity(newField, e))
1099  return SKIP;
1100  if (!this->requestLocality(&e, 1))
1101  return REQUEST;
1102  vertex = e;
1103  return OK;
1104  }
1105  virtual void apply()
1106  {
1107  /* the smoothing function used here is the average of the
1108  vertex value and neighboring vertex values, with the
1109  center vertex weighted equally as the neighbors */
1110  apf::Up edges;
1111  mesh->getUp(vertex, edges);
1112  apf::getComponents(field, vertex, 0, &sum[0]);
1113  for (int i = 0; i < edges.n; ++i)
1114  {
1115  apf::MeshEntity *ov = apf::getEdgeVertOppositeVert(
1116  mesh, edges.e[i], vertex);
1117  apf::getComponents(field, ov, 0, &value[0]);
1118  sum += value;
1119  }
1120  sum /= edges.n + 1;
1121  apf::setComponents(newField, vertex, 0, &sum[0]);
1122  ++nApplied;
1123  }
1124  apf::Field *field;
1125  apf::Field *newField;
1126  apf::MeshEntity *vertex;
1127  apf::DynamicVector sum;
1128  apf::DynamicVector value;
1129  apf::MeshTag *emptyTag;
1131 };
1132 
1133 static void SmoothField(apf::Field *f)
1134 {
1135  Smoother op(f);
1136  op.applyToDimension(0);
1137 }
1138 
1139 void getTargetError(apf::Mesh* m, apf::Field* errField, double &target_error,double totalError){
1140  //Implemented for 3D and for serial case only so far
1141  //Need to communicate target error in parallel
1142  logEvent("WARNING/ERROR:Parallel implementation is not completed yet",4);
1143  if(m->getDimension()==2){
1144  target_error = totalError/sqrt(m->count(m->getDimension()));
1145  if(PCU_Comm_Self()==0)
1146  std::cout<<"The estimated target error is "<<target_error<<std::endl;
1147  return;
1148  }
1149  apf::Field* interfaceField = m->findField("vof");
1150  apf::Field* targetField = apf::createField(m,"targetError",apf::SCALAR,apf::getVoronoiShape(m->getDimension(),1));
1151  apf::MeshEntity* ent;
1152  apf::MeshIterator* it = m->begin(m->getDimension());
1153  apf::MeshElement* element;
1154  apf::Element* vofElem;
1155  std::vector <double> errVect;
1156  while( (ent = m->iterate(it))){
1157  element = apf::createMeshElement(m, ent);
1158  vofElem = apf::createElement(interfaceField,element);
1159  double vofVal = apf::getScalar(vofElem,apf::Vector3(1./3.,1./3.,1./3.));
1160 
1161  if(vofVal < 0.9 && vofVal > 0.1){ //at the interface
1162  double errorValue = apf::getScalar(errField,ent,0);
1163  errVect.push_back(errorValue);
1164  apf::setScalar(targetField,ent,0,errorValue);
1165  }
1166  else{
1167  apf::setScalar(targetField,ent,0,0.0);
1168  }
1169  }
1170  m->end(it);
1171  if(errVect.size()==0){
1172  target_error = totalError/sqrt(m->count(m->getDimension()));
1173  }
1174  else{
1175  std::ofstream myfile;
1176  myfile.open("interfaceErrors.txt", std::ios::app );
1177  for(int i=0;i<errVect.size();i++){
1178  myfile << errVect[i]<<std::endl;
1179  }
1180  myfile.close();
1181  std::sort(errVect.begin(),errVect.end());
1182  int vectorSize = errVect.size();
1183  if(vectorSize %2 ==0){
1184  int idx1 = vectorSize/2-1;
1185  target_error = (errVect[idx1]+errVect[idx1+1])/2; //get average
1186  }
1187  else
1188  target_error = errVect[(vectorSize-1)/2];
1189  }
1190  if(PCU_Comm_Self()==0)
1191  std::cout<<"The estimated target error is "<<target_error<<std::endl;
1192  //std::abort();
1193 }
1194 
1196 //High level function that obtains the size scales and the size frames for ERM-based adapt and uses the computed total error
1197 {
1198 
1199  freeField(size_frame);
1200  freeField(size_scale);
1201  freeField(size_iso);
1202 
1203  //Initialize fields and needed types/variables
1204  apf::Field* errField;
1205  //apf::Mesh* m;
1206  if(hasERM)
1207  errField = m->findField("ErrorRegion");
1208  else if(hasVMS)
1209  errField = m->findField("VMSH1");
1210  assert(errField);
1211  //apf::Mesh *m = apf::getMesh(vmsErrH1);
1212  //apf::getMesh(errField);
1213  apf::MeshIterator *it;
1214  apf::MeshEntity *v;
1215  apf::MeshElement *element;
1216  apf::MeshEntity *reg;
1217  //size_iso = apf::createLagrangeField(m, "proteus_size", apf::SCALAR, 1);
1218  if(m->findField("errorSize"))
1219  apf::destroyField(m->findField("errorSize"));
1220  apf::Field *errorSize = apf::createLagrangeField(m, "errorSize", apf::SCALAR, 1);
1221 
1222  if (adapt_type_config == "anisotropic"){
1223  size_scale = apf::createLagrangeField(m, "proteus_size_scale", apf::VECTOR, 1);
1224  size_frame = apf::createLagrangeField(m, "proteus_size_frame", apf::MATRIX, 1);
1225  }
1226  apf::Field *errorSize_reg = apf::createField(m, "iso_size", apf::SCALAR, apf::getConstant(nsd));
1227  apf::Field *clipped_vtx = apf::createLagrangeField(m, "iso_clipped", apf::SCALAR, 1);
1228 
1229  //Get total number of elements
1230  int numel = 0;
1231  int nsd = m->getDimension();
1232  numel = m->count(nsd);
1233  PCU_Add_Ints(&numel, 1);
1234 
1235  //if target error is not specified, choose one based on interface or based on equidistribution assumption
1236  if(target_error==0){
1237  if(m->findField("vof")!=NULL)
1238  getTargetError(m,errField,target_error,err_total);
1239  else
1240  target_error = err_total/sqrt(m->count(nsd));
1241  }
1242 
1243  // Get domain volume
1244  // should only need to be computed once unless geometry is complex
1245  if (domainVolume == 0)
1246  {
1247  double volTotal = 0.0;
1248  it = m->begin(nsd);
1249  while (reg = m->iterate(it))
1250  {
1251  volTotal += apf::measure(m, reg);
1252  }
1253  PCU_Add_Doubles(&volTotal, 1);
1254  domainVolume = volTotal;
1255  assert(domainVolume > 0);
1256  }
1257  //compute the new size field over elements
1258  it = m->begin(nsd);
1259  double err_curr = 0.0;
1260  double errRho_curr = 0.0;
1261  double errRho_target = target_error / sqrt(domainVolume);
1262  apf::Vector3 err_vect;
1263  while (reg = m->iterate(it))
1264  {
1265  double h_old;
1266  double h_new;
1267  element = apf::createMeshElement(m, reg);
1268 
1269  if (m->getDimension() == 2)
1270  //h_old = sqrt(apf::measure(element) * 4 / sqrt(3));
1271  h_old = apf::computeLargestHeightInTri(m,reg);
1272  else
1273  //h_old = pow(apf::measure(element) * 6 * sqrt(2), 1.0 / 3.0); //edge of a regular tet
1274  h_old = apf::computeLargestHeightInTet(m,reg);
1275  //err_curr = apf::getScalar(vmsErrH1, reg, 0);
1276  err_curr = apf::getScalar(errField, reg, 0);
1277  //err_curr = err_vect[0];
1278  //errRho_curr = apf::getScalar(errRho_reg, reg, 0);
1279  //h_new = h_old*errRho_target/errRho_curr;
1280  //h_new = h_old*sqrt(apf::measure(element))/sqrt(domainVolume)*target_error/err_curr;
1281  //
1282  //error-to-size relationship should be different between anisotropic and isotropic cases
1283  //consider moving this to where size frames are computed to get aspect ratio info
1284  if (adapt_type_config == "anisotropic")
1285  if(target_error/err_curr <= 1)
1286  h_new = h_old * pow((target_error / err_curr),2.0/(2.0*(1.0)+1.0)); //refinement
1287  else
1288  h_new = h_old * pow((target_error / err_curr),2.0/(2.0*(1.0)+3.0)); //coarsening
1289  else //isotropic
1290  {
1291  if(target_error/err_curr <= 1)
1292  h_new = h_old * pow((target_error / err_curr),2.0/(2.0*(1.0)+nsd));
1293  else
1294  h_new = h_old * pow((target_error / err_curr),2.0/(2.0*(1.0)+nsd));
1295  //h_new = h_old * pow((target_error / err_curr),2.0/(2.0*(1.0)+nsd+1.0)); //extra 1.0 to slow down coarsening
1296  }
1297 
1298  apf::setScalar(errorSize_reg, reg, 0, h_new);
1299  apf::destroyMeshElement(element);
1300  }
1301  m->end(it);
1302 
1303  //Transfer size field from elements to vertices through averaging
1304  PCU_Comm_Begin();
1305  it = m->begin(0);
1306  while ((v = m->iterate(it)))
1307  {
1308  //averageToEntity(errorSize_reg, errorSize, v);
1309  //volumeAverageToEntity(errorSize_reg, errorSize, v);
1310  //errorAverageToEntity(errorSize_reg, errorSize,errField, v);
1311  minToEntity(errorSize_reg, errorSize, v);
1312  if(!m->isOwned(v))
1313  {
1314  apf::Copies remotes;
1315  m->getRemotes(v,remotes);
1316  int owningPart=m->getOwner(v);
1317  PCU_COMM_PACK(owningPart, remotes[owningPart]);
1318  double currentSize = apf::getScalar(errorSize,v,0);
1319  PCU_COMM_PACK(owningPart,currentSize);
1320  }
1321 
1322  }
1323  m->end(it);
1324  PCU_Comm_Send();
1325  while(PCU_Comm_Receive())
1326  {
1327  apf::MeshEntity* receivedEnt;
1328  double receivedSize;
1329  PCU_COMM_UNPACK(receivedEnt);
1330  PCU_COMM_UNPACK(receivedSize);
1331 
1332  double currentSize = apf::getScalar(errorSize,receivedEnt,0);
1333 
1334  if(receivedSize<currentSize)
1335  {
1336  apf::setScalar(errorSize,receivedEnt,0,receivedSize);
1337  }
1338 
1339  }
1340 
1341 
1342 
1343  //Get the anisotropic size frame
1344  if (adapt_type_config == "anisotropic")
1345  {
1346  if(comm_rank==0)
1347  std::cout<<"Entering anisotropic loop to compute size scales and frames\n";
1348  double eps_u = 0.002; //distance from the interface
1349 /*
1350  apf::Field *phif = m->findField("phi");
1351  apf::Field *gradphi = apf::recoverGradientByVolume(phif);
1352  apf::Field *grad2phi = apf::recoverGradientByVolume(gradphi);
1353 */
1354  apf::Field *speedF = extractSpeed(m->findField("velocity"));
1355  apf::Field *gradSpeed = apf::recoverGradientByVolume(speedF);
1356  apf::Field *grad2Speed = apf::recoverGradientByVolume(gradSpeed);
1357  //apf::Field *hess = computeHessianField(grad2phi);
1358  //apf::Field *curves = getCurves(hess, gradphi);
1359  //apf::Field* metricf = computeMetricField(gradphi,grad2phi,errorSize,eps_u);
1360  apf::Field *metricf = computeMetricField(gradSpeed, grad2Speed, errorSize, eps_u);
1361  apf::Field *frame_comps[3] = {apf::createLagrangeField(m, "frame_0", apf::VECTOR, 1), apf::createLagrangeField(m, "frame_1", apf::VECTOR, 1), apf::createLagrangeField(m, "frame_2", apf::VECTOR, 1)};
1362  //getERMSizeFrames(metricf, gradSpeed, frame_comps);
1363 
1364  //Set the size scale for vertices
1365  it = m->begin(0);
1366  apf::Vector3 scale;
1367  while ((v = m->iterate(it)))
1368  {
1369  double tempScale = apf::getScalar(errorSize, v, 0);
1370  if (tempScale < hmin)
1371  apf::setScalar(clipped_vtx, v, 0, -1);
1372  else if (tempScale > hmax)
1373  apf::setScalar(clipped_vtx, v, 0, 1);
1374  else
1375  apf::setScalar(clipped_vtx, v, 0, 0);
1376  clamp(tempScale, hmin, hmax);
1377  apf::setScalar(errorSize,v,0,tempScale);
1378  }
1379  it = m->begin(0);
1380  while( (v = m->iterate(it)) ){
1381  double phi;// = apf::getScalar(phif, v, 0);
1382  apf::Vector3 curve;
1383  //apf::getVector(curves, v, 0, curve);
1384 
1385  //metricf is the hessian
1386  apf::Matrix3x3 metric;
1387  apf::getMatrix(metricf, v, 0, metric);
1388 
1389  apf::Vector3 eigenVectors[3];
1390  double eigenValues[3];
1391  apf::eigen(metric, eigenVectors, eigenValues);
1392  // Sort the eigenvalues and corresponding vectors
1393  // Larger eigenvalues means a need for a finer mesh
1394  SortingStruct ssa[3];
1395  for (int i = 0; i < 3; ++i)
1396  {
1397  ssa[i].v = eigenVectors[i];
1398  ssa[i].wm = std::fabs(eigenValues[i]);
1399  }
1400  std::sort(ssa, ssa + 3);
1401 
1402  assert(ssa[2].wm >= ssa[1].wm);
1403  assert(ssa[1].wm >= ssa[0].wm);
1404 
1405  double lambda[3] = {ssa[2].wm, ssa[1].wm, ssa[0].wm};
1406 
1407  scaleFormulaERM(phi, hmin, hmax, apf::getScalar(errorSize, v, 0), curve, lambda, eps_u, scale,nsd,maxAspect);
1408  apf::setVector(size_scale, v, 0, scale);
1409  //get frames
1410 
1411  apf::Matrix3x3 frame(1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0);
1412 
1413  //get eigen values and eigenvectors from hessian
1414  double firstEigenvalue = ssa[2].wm;
1415  assert(firstEigenvalue > 1e-12);
1416  frame[0] = ssa[2].v;
1417  frame[1] = ssa[1].v;
1418  frame[2] = ssa[0].v;
1419 
1420  //normalize eigenvectors
1421  for (int i = 0; i < 3; ++i)
1422  frame[i] = frame[i].normalize();
1423  frame = apf::transpose(frame);
1424  apf::setMatrix(size_frame, v, 0, frame);
1425 
1426  }
1427  m->end(it);
1428 
1429  //Do simple size and aspect ratio grading
1431  if(comm_rank==0)
1432  std::cout<<"Finished grading size 0\n";
1434  if(comm_rank==0)
1435  std::cout<<"Finished grading size 1\n";
1437  if(comm_rank==0)
1438  std::cout<<"Finished grading size 2\n";
1439 
1440  apf::synchronize(size_scale);
1441 
1442  //apf::destroyField(gradphi);
1443  //apf::destroyField(grad2phi);
1444  //apf::destroyField(curves);
1445  //apf::destroyField(hess);
1446 
1447  if(logging_config=="on"){
1448  char namebuffer[20];
1449  sprintf(namebuffer,"pumi_preadapt_aniso_%i",nAdapt);
1450  apf::writeVtkFiles(namebuffer, m);
1451  }
1452 
1453  apf::destroyField(metricf);
1454  apf::destroyField(frame_comps[0]);
1455  apf::destroyField(frame_comps[1]);
1456  apf::destroyField(frame_comps[2]);
1457  apf::destroyField(speedF);
1458  apf::destroyField(gradSpeed);
1459  apf::destroyField(grad2Speed);
1460  }
1461  else
1462  {
1463  it = m->begin(0);
1464  while ((v = m->iterate(it)))
1465  {
1466  double tempScale = apf::getScalar(errorSize, v, 0);
1467  if (tempScale < hmin)
1468  apf::setScalar(clipped_vtx, v, 0, -1);
1469  else if (tempScale > hmax)
1470  apf::setScalar(clipped_vtx, v, 0, 1);
1471  else
1472  apf::setScalar(clipped_vtx, v, 0, 0);
1473  clamp(tempScale, hmin, hmax);
1474  apf::setScalar(errorSize, v, 0, tempScale);
1475  }
1476  //gradeMesh();
1477  apf::synchronize(errorSize);
1478  m->end(it);
1479  if (target_element_count != 0)
1480  {
1481  sam::scaleIsoSizeField(errorSize, target_element_count);
1482  clampField(errorSize, hmin, hmax);
1483  //gradeMesh();
1484  //SmoothField(errorSize);
1485  }
1486  sizeFieldList.push(errorSize);
1487  }
1488 
1489  //Destroy locally required fields
1490  apf::destroyField(errorSize_reg);
1491  apf::destroyField(clipped_vtx);
1492  if(comm_rank==0)
1493  std::cout<<"Finished Size Field\n";
1494  return 0;
1495 }
1496 
1498 //Function that tests MeshAdapt by generating an isotropic sizefield based on hmin
1499 {
1500  freeField(size_iso);
1501  size_iso = apf::createLagrangeField(m, "proteus_size",apf::SCALAR,1);
1502  apf::MeshIterator* it = m->begin(0);
1503  apf::MeshEntity* v;
1504  while(v = m->iterate(it)){
1505  double phi = hmin;
1506  clamp(phi,hmin,hmax);
1507  apf::setScalar(size_iso,v,0,phi);
1508  }
1509  return 0;
1510 }
1511 
1512 int gradeSizeModify(apf::Mesh* m, double gradingFactor,
1513  double size[2], apf::Adjacent edgAdjVert,
1514  apf::Adjacent vertAdjEdg,
1515  std::queue<apf::MeshEntity*> &markedEdges,
1516  apf::MeshTag* isMarked,
1517  int fieldType,
1518  int vecPos, //which idx of sizeVec to modify
1519  int idxFlag)
1520 
1521 //General function to actually modify sizes
1522 {
1523  //Determine a switching scheme depending on which vertex needs a modification
1524  int idx1,idx2;
1525  if(idxFlag == 0){
1526  idx1=0;
1527  idx2=1;
1528  }
1529  else{
1530  idx1=1;
1531  idx2 = 0;
1532  }
1533 
1534  int marker[3] = {0,1,0};
1535  double marginVal = 0.01;
1536  int needsParallel=0;
1537 
1538  if(fieldType == apf::SCALAR){
1539  apf::Field* size_iso = m->findField("proteus_size");
1540 
1541  if(size[idx1]>(gradingFactor*size[idx2])*(1+marginVal))
1542  {
1543  if(m->isOwned(edgAdjVert[idx1]))
1544  {
1545  size[idx1] = gradingFactor*size[idx2];
1546  apf::setScalar(size_iso,edgAdjVert[idx1],0,size[idx1]);
1547  m->getAdjacent(edgAdjVert[idx1], 1, vertAdjEdg);
1548  for (std::size_t i=0; i<vertAdjEdg.getSize();++i){
1549  m->getIntTag(vertAdjEdg[i],isMarked,&marker[2]);
1550  //if edge is not already marked
1551  if(!marker[2]){
1552  m->setIntTag(vertAdjEdg[i],isMarked,&marker[1]);
1553  markedEdges.push(vertAdjEdg[i]);
1554  }
1555  }
1556  } //end isOwned
1557  else
1558  { //Pack information to owning processor
1559  needsParallel=1;
1560  apf::Copies remotes;
1561  m->getRemotes(edgAdjVert[idx1],remotes);
1562  double newSize = gradingFactor*size[idx2];
1563  int owningPart=m->getOwner(edgAdjVert[idx1]);
1564  PCU_COMM_PACK(owningPart, remotes[owningPart]);
1565  PCU_COMM_PACK(owningPart,newSize);
1566  }
1567  }
1568 
1569  }//end if apf::SCALAR
1570  else{
1571  apf::Field* size_scale = m->findField("proteus_size_scale");
1572  apf::Vector3 sizeVec;
1573  if(size[idx1]>(gradingFactor*size[idx2])*(1+marginVal)){
1574  size[idx1] = gradingFactor*size[idx2];
1575  apf::getVector(size_scale,edgAdjVert[idx1],0,sizeVec);
1576  if(vecPos > 0){
1577  sizeVec[vecPos] = size[idx1]*sizeVec[0]; //realize the new aspect ratio
1578  }
1579  else{
1580  sizeVec[0] = size[idx1];
1581  }
1582  apf::setVector(size_scale,edgAdjVert[idx1],0,sizeVec);
1583  m->getAdjacent(edgAdjVert[idx1], 1, vertAdjEdg);
1584  for (std::size_t i=0; i<vertAdjEdg.getSize();++i){
1585  m->getIntTag(vertAdjEdg[i],isMarked,&marker[2]);
1586  //if edge is not already marked
1587  if(!marker[2]){
1588  m->setIntTag(vertAdjEdg[i],isMarked,&marker[1]);
1589  markedEdges.push(vertAdjEdg[i]);
1590  }
1591  }
1592  }
1593  }
1594  return needsParallel;
1595 }
1596 
1597 void markEdgesInitial(apf::Mesh* m, std::queue<apf::MeshEntity*> &markedEdges,double gradingFactor)
1598 //Function used to initially determine which edges need to be considered for gradation
1599 {
1600  //marker structure for 0) not marked 1) marked 2)storage
1601  int marker[3] = {0,1,0};
1602 
1603  double size[2];
1604  apf::MeshTag* isMarked = m->findTag("isMarked");
1605  apf::Field* size_iso = m->findField("proteus_size");
1606  apf::Adjacent edgAdjVert;
1607  apf::MeshEntity* edge;
1608  apf::MeshIterator* it = m->begin(1);
1609  while((edge=m->iterate(it))){
1610  m->getAdjacent(edge, 0, edgAdjVert);
1611  for (std::size_t i=0; i < edgAdjVert.getSize(); ++i){
1612  size[i]=apf::getScalar(size_iso,edgAdjVert[i],0);
1613  }
1614  if( (size[0] > gradingFactor*size[1]) || (size[1] > gradingFactor*size[0]) ){
1615  //add edge to a queue
1616  markedEdges.push(edge);
1617  //tag edge to indicate that it is part of queue
1618  m->setIntTag(edge,isMarked,&marker[1]);
1619  }
1620  else{
1621  m->setIntTag(edge,isMarked,&marker[0]);
1622  }
1623  }
1624  m->end(it);
1625 }
1626 
1627 int serialGradation(apf::Mesh* m, std::queue<apf::MeshEntity*> &markedEdges,double gradingFactor)
1628 //Function used loop over the mesh edge queue for gradation and modify the sizes
1629 {
1630  double size[2];
1631  //marker structure for 0) not marked 1) marked 2)storage
1632  int marker[3] = {0,1,0};
1633  apf::MeshTag* isMarked = m->findTag("isMarked");
1634  apf::Field* size_iso = m->findField("proteus_size");
1635  apf::Adjacent edgAdjVert;
1636  apf::Adjacent vertAdjEdg;
1637  apf::MeshEntity* edge;
1638  apf::MeshIterator* it = m->begin(1);
1639  int needsParallel=0;
1640 
1641  //perform serial gradation while packing necessary info for parallel
1642  while(!markedEdges.empty()){
1643  edge = markedEdges.front();
1644  m->getAdjacent(edge, 0, edgAdjVert);
1645  for (std::size_t i=0; i < edgAdjVert.getSize(); ++i){
1646  size[i] = apf::getScalar(size_iso,edgAdjVert[i],0);
1647  }
1648 
1649  needsParallel+=gradeSizeModify(m, gradingFactor, size, edgAdjVert,
1650  vertAdjEdg, markedEdges, isMarked, apf::SCALAR,0, 0);
1651  needsParallel+=gradeSizeModify(m, gradingFactor, size, edgAdjVert,
1652  vertAdjEdg, markedEdges, isMarked, apf::SCALAR,0, 1);
1653 
1654  m->setIntTag(edge,isMarked,&marker[0]);
1655  markedEdges.pop();
1656  }
1657  return needsParallel;
1658 }
1659 
1660 int MeshAdaptPUMIDrvr::gradeMesh(double gradationFactor)
1661 //Function to grade isotropic mesh through comparison of edge vertex size ratios
1662 //This implementation accounts for parallel meshes as well
1663 //First do serial gradation.
1664 //If a shared entity has its size modified, then send new size to owning copy.
1665 //After full loop over entities, have owning copy take minimum of all sizes received
1666 //Flag adjacent entities to owning copy.
1667 //Communicate to remote copies that a size was modified, and flag adjacent edges to remote copies for further gradation
1668 {
1669  //
1670  logEvent("Starting grading",4);
1671  apf::MeshEntity* edge;
1672  apf::Adjacent edgAdjVert;
1673  apf::Adjacent vertAdjEdg;
1674  double size[2];
1675  std::queue<apf::MeshEntity*> markedEdges;
1676  apf::MeshTag* isMarked = m->createIntTag("isMarked",1);
1677 
1678  //marker structure for 0) not marked 1) marked 2)storage
1679  int marker[3] = {0,1,0};
1680 
1681  apf::MeshIterator* it;
1682  markEdgesInitial(m,markedEdges,gradingFactor);
1683 
1684  int needsParallel=1;
1685  int nCount=1;
1686  while(needsParallel)
1687  {
1688  PCU_Comm_Begin();
1689  needsParallel = serialGradation(m,markedEdges,gradingFactor);
1690 
1691  PCU_Add_Ints(&needsParallel,1);
1692  if(comm_rank==0)
1693  std::cerr<<"Sending size info for gradation"<<std::endl;
1694  PCU_Comm_Send();
1695 
1696  apf::MeshEntity* ent;
1697  double receivedSize;
1698  double currentSize;
1699  double newSize;
1700 
1701  //Need a container to get all entitites that need to be updated on remotes
1702  std::queue<apf::MeshEntity*> updateRemoteVertices;
1703 
1704  apf::Copies remotes;
1705  //owning copies are receiving
1706  while(PCU_Comm_Receive())
1707  {
1708  PCU_COMM_UNPACK(ent);
1709  PCU_COMM_UNPACK(receivedSize);
1710 
1711  if(!m->isOwned(ent)){
1712  logEvent("THERE WAS AN ERROR",4);
1713  std::exit(1);
1714  }
1715 
1716  currentSize = apf::getScalar(size_iso,ent,0);
1717  newSize = std::min(receivedSize,currentSize);
1718  apf::setScalar(size_iso,ent,0,newSize);
1719 
1720  //add adjacent edges into Q
1721  m->getAdjacent(ent, 1, vertAdjEdg);
1722  for (std::size_t i=0; i<vertAdjEdg.getSize();++i)
1723  {
1724  edge = vertAdjEdg[i];
1725  m->getIntTag(vertAdjEdg[i],isMarked,&marker[2]);
1726  if(!marker[2])
1727  {
1728  markedEdges.push(edge);
1729  //tag edge to indicate that it is part of queue
1730  m->setIntTag(edge,isMarked,&marker[1]);
1731  }
1732  }
1733  updateRemoteVertices.push(ent);
1734  }
1735 
1736  PCU_Comm_Begin();
1737 
1738  while(!updateRemoteVertices.empty())
1739  {
1740  ent = updateRemoteVertices.front();
1741  //get remote copies and send updated mesh sizes
1742  m->getRemotes(ent,remotes);
1743  currentSize = apf::getScalar(size_iso,ent,0);
1744  for(apf::Copies::iterator iter=remotes.begin(); iter!=remotes.end();++iter)
1745  {
1746  PCU_COMM_PACK(iter->first, iter->second);
1747  }
1748  updateRemoteVertices.pop();
1749  }
1750 
1751  PCU_Comm_Send();
1752  //while remote copies are receiving
1753  while(PCU_Comm_Receive())
1754  {
1755  //unpack
1756  PCU_COMM_UNPACK(ent);
1757  //PCU_COMM_UNPACK(receivedSize);
1758  assert(!m->isOwned(ent));
1759 
1760  if(m->isOwned(ent)){
1761  logEvent("Problem occurred",4);
1762  std::exit(1);
1763  }
1764 
1765  //add adjacent edges into Q
1766  m->getAdjacent(ent, 1, vertAdjEdg);
1767  for (std::size_t i=0; i<vertAdjEdg.getSize();++i)
1768  {
1769  edge = vertAdjEdg[i];
1770  m->getIntTag(vertAdjEdg[i],isMarked,&marker[2]);
1771  if(!marker[2])
1772  {
1773  markedEdges.push(edge);
1774  //tag edge to indicate that it is part of queue
1775  m->setIntTag(edge,isMarked,&marker[1]);
1776  }
1777  }
1778  }
1779  apf::synchronize(size_iso);
1780 
1781  } //end outer while
1782 
1783  //Cleanup of edge marker field
1784  it = m->begin(1);
1785  while((edge=m->iterate(it))){
1786  m->removeTag(edge,isMarked);
1787  }
1788  m->end(it);
1789  m->destroyTag(isMarked);
1790 
1791  //apf::synchronize(size_iso);
1792  logEvent("Completed grading",4);
1793  return needsParallel;
1794 }
1795 
1796 void gradeAnisoMesh(apf::Mesh* m)
1797 //Function to grade anisotropic mesh through comparison of edge vertex aspect ratios and minimum sizes
1798 //For simplicity, we do not bother with accounting for entities across partitions
1799 {
1800  //
1801  //if(comm_rank==0)
1802  // std::cout<<"Starting anisotropic grading\n";
1803  apf::MeshIterator* it = m->begin(1);
1804  apf::MeshEntity* edge;
1805  apf::Adjacent edgAdjVert;
1806  apf::Adjacent vertAdjEdg;
1807  double gradingFactor = 1.3;
1808  double size[2];
1809  apf::Vector3 sizeVec;
1810  std::queue<apf::MeshEntity*> markedEdges;
1811  apf::MeshTag* isMarked = m->createIntTag("isMarked",1);
1812  apf::Field* size_scale = m->findField("proteus_size_scale");
1813 
1814  //marker structure for 0) not marked 1) marked 2)storage
1815  int marker[3] = {0,1,0};
1816 
1817  while((edge=m->iterate(it))){
1818  m->getAdjacent(edge, 0, edgAdjVert);
1819  for (std::size_t i=0; i < edgAdjVert.getSize(); ++i){
1820  apf::getVector(size_scale,edgAdjVert[i],0,sizeVec);
1821  size[i]=sizeVec[0];
1822  }
1823  if( (size[0] > gradingFactor*size[1]) || (size[1] > gradingFactor*size[0]) ){
1824  //add edge to a queue
1825  markedEdges.push(edge);
1826  //tag edge to indicate that it is part of queue
1827  m->setIntTag(edge,isMarked,&marker[1]);
1828 
1829  }
1830  else{
1831  m->setIntTag(edge,isMarked,&marker[0]);
1832  }
1833  }
1834  m->end(it);
1835  while(!markedEdges.empty()){
1836  edge = markedEdges.front();
1837  m->getAdjacent(edge, 0, edgAdjVert);
1838  for (std::size_t i=0; i < edgAdjVert.getSize(); ++i){
1839  apf::getVector(size_scale,edgAdjVert[i],0,sizeVec);
1840  size[i]=sizeVec[0];
1841  }
1842  gradeSizeModify(m, gradingFactor, size, edgAdjVert,
1843  vertAdjEdg, markedEdges, isMarked, apf::VECTOR,0, 0);
1844  gradeSizeModify(m, gradingFactor, size, edgAdjVert,
1845  vertAdjEdg, markedEdges, isMarked, apf::VECTOR,0, 1);
1846 
1847 /*
1848  if(size[0]>gradingFactor*size[1]){
1849  size[0] = gradingFactor*size[1];
1850  apf::getVector(size_scale,edgAdjVert[0],0,sizeVec);
1851  sizeVec[0] = size[0];
1852  apf::setVector(size_scale,edgAdjVert[0],0,sizeVec);
1853  m->getAdjacent(edgAdjVert[0], 1, vertAdjEdg);
1854  for (std::size_t i=0; i<vertAdjEdg.getSize();++i){
1855  m->getIntTag(vertAdjEdg[i],isMarked,&marker[2]);
1856  //if edge is not already marked
1857  if(!marker[2]){
1858  m->setIntTag(vertAdjEdg[i],isMarked,&marker[1]);
1859  markedEdges.push(vertAdjEdg[i]);
1860  }
1861  }
1862  }
1863  if(size[1]>gradingFactor*size[0]){
1864  size[1] = gradingFactor*size[0];
1865  apf::getVector(size_scale,edgAdjVert[1],0,sizeVec);
1866  sizeVec[0] = size[1];
1867  apf::setVector(size_scale,edgAdjVert[1],0,sizeVec);
1868  m->getAdjacent(edgAdjVert[1], 1, vertAdjEdg);
1869  for (std::size_t i=0; i<vertAdjEdg.getSize();++i){
1870  m->getIntTag(vertAdjEdg[i],isMarked,&marker[2]);
1871  //if edge is not already marked
1872  if(!marker[2]){
1873  m->setIntTag(vertAdjEdg[i],isMarked,&marker[1]);
1874  markedEdges.push(vertAdjEdg[i]);
1875  }
1876  }
1877  }
1878 */
1879  m->setIntTag(edge,isMarked,&marker[0]);
1880  markedEdges.pop();
1881  }
1882  it = m->begin(1);
1883  while((edge=m->iterate(it))){
1884  m->removeTag(edge,isMarked);
1885  }
1886  m->end(it);
1887  m->destroyTag(isMarked);
1888  apf::synchronize(size_scale);
1889  //if(comm_rank==0)
1890  // std::cout<<"Completed minimum size grading\n";
1891 }
1892 
1893 void gradeAspectRatio(apf::Mesh* m,int idx)
1894 //Function to grade anisotropic mesh through comparison of edge vertex aspect ratios and minimum sizes
1895 //For simplicity, we do not bother with accounting for entities across partitions
1896 {
1897  std::cout<<"Entered function\n";
1898  apf::MeshIterator* it = m->begin(1);
1899  apf::MeshEntity* edge;
1900  apf::Adjacent edgAdjVert;
1901  apf::Adjacent vertAdjEdg;
1902  double gradingFactor = 1.3;
1903  double size[2];
1904  apf::Vector3 sizeVec;
1905  std::queue<apf::MeshEntity*> markedEdges;
1906  apf::MeshTag* isMarked = m->createIntTag("isMarked",1);
1907  apf::Field* size_scale = m->findField("proteus_size_scale");
1908 
1909  //marker structure for 0) not marked 1) marked 2)storage
1910  int marker[3] = {0,1,0};
1911 
1912  while((edge=m->iterate(it))){
1913  m->getAdjacent(edge, 0, edgAdjVert);
1914  for (std::size_t i=0; i < edgAdjVert.getSize(); ++i){
1915  apf::getVector(size_scale,edgAdjVert[i],0,sizeVec);
1916  size[i]=sizeVec[idx]/sizeVec[0];
1917  }
1918  if( (size[0] > gradingFactor*size[1]) || (size[1] > gradingFactor*size[0]) ){
1919  //add edge to a queue
1920  markedEdges.push(edge);
1921  //tag edge to indicate that it is part of queue
1922  m->setIntTag(edge,isMarked,&marker[1]);
1923 
1924  }
1925  else{
1926  m->setIntTag(edge,isMarked,&marker[0]);
1927  }
1928  }
1929  m->end(it);
1930 
1931  std::cout<<"Got queue of size "<<markedEdges.size()<<std::endl;
1932  while(!markedEdges.empty()){
1933  edge = markedEdges.front();
1934  m->getAdjacent(edge, 0, edgAdjVert);
1935  for (std::size_t i=0; i < edgAdjVert.getSize(); ++i){
1936  apf::getVector(size_scale,edgAdjVert[i],0,sizeVec);
1937  size[i]=sizeVec[idx]/sizeVec[0];
1938  }
1939  gradeSizeModify(m, gradingFactor, size, edgAdjVert,
1940  vertAdjEdg, markedEdges, isMarked, apf::VECTOR, idx, 0);
1941  gradeSizeModify(m, gradingFactor, size, edgAdjVert,
1942  vertAdjEdg, markedEdges, isMarked, apf::VECTOR, idx, 1);
1943 
1944  m->setIntTag(edge,isMarked,&marker[0]);
1945  markedEdges.pop();
1946  }
1947  it = m->begin(1);
1948  while((edge=m->iterate(it))){
1949  m->removeTag(edge,isMarked);
1950  }
1951  m->end(it);
1952  m->destroyTag(isMarked);
1953  apf::synchronize(size_scale);
1954 }
1955 
1956 void gradeAnisoMesh(apf::Mesh* m,double gradingFactor)
1957 //Function to grade anisotropic mesh through comparison of edge vertex aspect ratios and minimum sizes
1958 //For simplicity, we do not bother with accounting for entities across partitions
1959 {
1960  //
1961  //if(comm_rank==0)
1962  // std::cout<<"Starting anisotropic grading\n";
1963  apf::MeshIterator* it = m->begin(1);
1964  apf::MeshEntity* edge;
1965  apf::Adjacent edgAdjVert;
1966  apf::Adjacent vertAdjEdg;
1967  //double gradingFactor = 1.3;
1968  double size[2];
1969  apf::Vector3 sizeVec;
1970  std::queue<apf::MeshEntity*> markedEdges;
1971  apf::MeshTag* isMarked = m->createIntTag("isMarked",1);
1972  apf::Field* size_scale = m->findField("proteus_size_scale");
1973 
1974  //marker structure for 0) not marked 1) marked 2)storage
1975  int marker[3] = {0,1,0};
1976 
1977  while((edge=m->iterate(it))){
1978  m->getAdjacent(edge, 0, edgAdjVert);
1979  for (std::size_t i=0; i < edgAdjVert.getSize(); ++i){
1980  apf::getVector(size_scale,edgAdjVert[i],0,sizeVec);
1981  size[i]=sizeVec[0];
1982  }
1983  if( (size[0] > gradingFactor*size[1]) || (size[1] > gradingFactor*size[0]) ){
1984  //add edge to a queue
1985  markedEdges.push(edge);
1986  //tag edge to indicate that it is part of queue
1987  m->setIntTag(edge,isMarked,&marker[1]);
1988 
1989  }
1990  else{
1991  m->setIntTag(edge,isMarked,&marker[0]);
1992  }
1993  }
1994  m->end(it);
1995  while(!markedEdges.empty()){
1996  edge = markedEdges.front();
1997  m->getAdjacent(edge, 0, edgAdjVert);
1998  for (std::size_t i=0; i < edgAdjVert.getSize(); ++i){
1999  apf::getVector(size_scale,edgAdjVert[i],0,sizeVec);
2000  size[i]=sizeVec[0];
2001  }
2002  gradeSizeModify(m, gradingFactor, size, edgAdjVert,
2003  vertAdjEdg, markedEdges, isMarked, apf::VECTOR,0, 0);
2004  gradeSizeModify(m, gradingFactor, size, edgAdjVert,
2005  vertAdjEdg, markedEdges, isMarked, apf::VECTOR,0, 1);
2006 
2007 /*
2008  if(size[0]>gradingFactor*size[1]){
2009  size[0] = gradingFactor*size[1];
2010  apf::getVector(size_scale,edgAdjVert[0],0,sizeVec);
2011  sizeVec[0] = size[0];
2012  apf::setVector(size_scale,edgAdjVert[0],0,sizeVec);
2013  m->getAdjacent(edgAdjVert[0], 1, vertAdjEdg);
2014  for (std::size_t i=0; i<vertAdjEdg.getSize();++i){
2015  m->getIntTag(vertAdjEdg[i],isMarked,&marker[2]);
2016  //if edge is not already marked
2017  if(!marker[2]){
2018  m->setIntTag(vertAdjEdg[i],isMarked,&marker[1]);
2019  markedEdges.push(vertAdjEdg[i]);
2020  }
2021  }
2022  }
2023  if(size[1]>gradingFactor*size[0]){
2024  size[1] = gradingFactor*size[0];
2025  apf::getVector(size_scale,edgAdjVert[1],0,sizeVec);
2026  sizeVec[0] = size[1];
2027  apf::setVector(size_scale,edgAdjVert[1],0,sizeVec);
2028  m->getAdjacent(edgAdjVert[1], 1, vertAdjEdg);
2029  for (std::size_t i=0; i<vertAdjEdg.getSize();++i){
2030  m->getIntTag(vertAdjEdg[i],isMarked,&marker[2]);
2031  //if edge is not already marked
2032  if(!marker[2]){
2033  m->setIntTag(vertAdjEdg[i],isMarked,&marker[1]);
2034  markedEdges.push(vertAdjEdg[i]);
2035  }
2036  }
2037  }
2038 */
2039  m->setIntTag(edge,isMarked,&marker[0]);
2040  markedEdges.pop();
2041  }
2042  it = m->begin(1);
2043  while((edge=m->iterate(it))){
2044  m->removeTag(edge,isMarked);
2045  }
2046  m->end(it);
2047  m->destroyTag(isMarked);
2048  apf::synchronize(size_scale);
2049  //if(comm_rank==0)
2050  // std::cout<<"Completed minimum size grading\n";
2051 }
2052 
2053 void gradeAspectRatio(apf::Mesh* m,int idx,double gradingFactor)
2054 //Function to grade anisotropic mesh through comparison of edge vertex aspect ratios and minimum sizes
2055 //For simplicity, we do not bother with accounting for entities across partitions
2056 {
2057  if(PCU_Comm_Self()==0)
2058  std::cout<<"Entered function\n";
2059  apf::MeshIterator* it = m->begin(1);
2060  apf::MeshEntity* edge;
2061  apf::Adjacent edgAdjVert;
2062  apf::Adjacent vertAdjEdg;
2063  //double gradingFactor = 1.3;
2064  double size[2];
2065  apf::Vector3 sizeVec;
2066  std::queue<apf::MeshEntity*> markedEdges;
2067  apf::MeshTag* isMarked = m->createIntTag("isMarked",1);
2068  apf::Field* size_scale = m->findField("proteus_size_scale");
2069 
2070  //marker structure for 0) not marked 1) marked 2)storage
2071  int marker[3] = {0,1,0};
2072 
2073  while((edge=m->iterate(it))){
2074  m->getAdjacent(edge, 0, edgAdjVert);
2075  for (std::size_t i=0; i < edgAdjVert.getSize(); ++i){
2076  apf::getVector(size_scale,edgAdjVert[i],0,sizeVec);
2077  size[i]=sizeVec[idx]/sizeVec[0];
2078  }
2079  if( (size[0] > gradingFactor*size[1]) || (size[1] > gradingFactor*size[0]) ){
2080  //add edge to a queue
2081  markedEdges.push(edge);
2082  //tag edge to indicate that it is part of queue
2083  m->setIntTag(edge,isMarked,&marker[1]);
2084 
2085  }
2086  else{
2087  m->setIntTag(edge,isMarked,&marker[0]);
2088  }
2089  }
2090  m->end(it);
2091 
2092  if(PCU_Comm_Self()==0)
2093  std::cout<<"Got queue of size "<<markedEdges.size()<<std::endl;
2094  while(!markedEdges.empty()){
2095  edge = markedEdges.front();
2096  m->getAdjacent(edge, 0, edgAdjVert);
2097  for (std::size_t i=0; i < edgAdjVert.getSize(); ++i){
2098  apf::getVector(size_scale,edgAdjVert[i],0,sizeVec);
2099  size[i]=sizeVec[idx]/sizeVec[0];
2100  }
2101  gradeSizeModify(m, gradingFactor, size, edgAdjVert,
2102  vertAdjEdg, markedEdges, isMarked, apf::VECTOR, idx, 0);
2103  gradeSizeModify(m, gradingFactor, size, edgAdjVert,
2104  vertAdjEdg, markedEdges, isMarked, apf::VECTOR, idx, 1);
2105 
2106  m->setIntTag(edge,isMarked,&marker[0]);
2107  markedEdges.pop();
2108  }
2109  it = m->begin(1);
2110  while((edge=m->iterate(it))){
2111  m->removeTag(edge,isMarked);
2112  }
2113  m->end(it);
2114  m->destroyTag(isMarked);
2115  apf::synchronize(size_scale);
2116 }
MeshAdaptPUMIDrvr::hasERM
bool hasERM
Definition: MeshAdaptPUMI.h:105
MeshAdaptPUMIDrvr::hPhi
double hPhi
Definition: MeshAdaptPUMI.h:91
Smoother::setEntity
virtual Outcome setEntity(apf::MeshEntity *e)
Definition: SizeField.cpp:1096
Smoother
Definition: SizeField.cpp:1081
MeshAdaptPUMIDrvr::hmin
double hmin
Definition: MeshAdaptPUMI.h:91
MeshAdaptPUMI.h
gradeAnisoMesh
void gradeAnisoMesh(apf::Mesh *m, double gradingFactor)
Definition: SizeField.cpp:1956
MeshAdaptPUMIDrvr::predictiveInterfacePropagation
void predictiveInterfacePropagation()
Definition: SizeField.cpp:297
MeshAdaptPUMIDrvr::setSphereSizeField
int setSphereSizeField()
Definition: SizeField.cpp:61
MeshAdaptPUMIDrvr::nsd
int nsd
Definition: MeshAdaptPUMI.h:96
Smoother::value
apf::DynamicVector value
Definition: SizeField.cpp:1128
MeshAdaptPUMIDrvr::hmax
double hmax
Definition: MeshAdaptPUMI.h:91
f
Double f
Definition: Headers.h:64
MeshAdaptPUMIDrvr::calculateAnisoSizeField
int calculateAnisoSizeField()
Definition: SizeField.cpp:1055
errorAverageToEntity
void errorAverageToEntity(apf::Field *ef, apf::Field *vf, apf::Field *err, apf::MeshEntity *ent)
Definition: SizeField.cpp:618
edgeWalkerInfo::trackerTag
apf::MeshTag * trackerTag
Definition: SizeField.cpp:185
Smoother::vertex
apf::MeshEntity * vertex
Definition: SizeField.cpp:1126
getTargetError
void getTargetError(apf::Mesh *m, apf::Field *errField, double &target_error, double totalError)
Definition: SizeField.cpp:1139
s
Double s
Definition: Headers.h:84
gradeAspectRatio
void gradeAspectRatio(apf::Mesh *m, int idx, double gradingFactor)
Definition: SizeField.cpp:2053
MeshAdaptPUMIDrvr::localNumber
int localNumber(apf::MeshEntity *e)
Definition: MeshConverter.cpp:105
SortingStruct::v
apf::Vector3 v
Definition: SizeField.cpp:925
Smoother::Smoother
Smoother(apf::Field *f)
Definition: SizeField.cpp:1082
MeshAdaptPUMIDrvr::getERMSizeField
int getERMSizeField(double err_total)
Definition: SizeField.cpp:1195
apf
Definition: MeshConverter.cpp:639
edgeWalkerInfo
Definition: SizeField.cpp:180
phi
Double phi
Definition: Headers.h:76
serialGradation
int serialGradation(apf::Mesh *m, std::queue< apf::MeshEntity * > &markedEdges, double gradingFactor)
Definition: SizeField.cpp:1627
MeshAdaptPUMIDrvr::hasVMS
bool hasVMS
Definition: MeshAdaptPUMI.h:104
MeshAdaptPUMIDrvr::maxAspect
int maxAspect
Definition: MeshAdaptPUMI.h:97
edgeWalkerInfo::direction
double direction
Definition: SizeField.cpp:183
edgeWalkerInfo::edgeID
int edgeID
Definition: SizeField.cpp:187
Smoother::~Smoother
~Smoother()
Definition: SizeField.cpp:1091
MeshAdaptPUMIDrvr::N_interface_band
double N_interface_band
Definition: MeshAdaptPUMI.h:100
min
#define min(a, b)
Definition: jf.h:71
MeshAdaptPUMIDrvr::gradingFactor
double gradingFactor
Definition: MeshAdaptPUMI.h:101
v
Double v
Definition: Headers.h:95
Smoother::newField
apf::Field * newField
Definition: SizeField.cpp:1125
markEdgesInitial
void markEdgesInitial(apf::Mesh *m, std::queue< apf::MeshEntity * > &markedEdges, double gradingFactor)
Definition: SizeField.cpp:1597
intersectsInterface
int intersectsInterface(apf::MeshEntity *edge, apf::Field *levelSet)
Definition: SizeField.cpp:164
SortingStruct::wm
double wm
Definition: SizeField.cpp:926
edgeLength
double edgeLength(int nL, int nR, const double *nodeArray)
Definition: mesh.cpp:3484
c
Double c
Definition: Headers.h:54
SortingStruct
Definition: SizeField.cpp:924
gradeSizeModify
int gradeSizeModify(apf::Mesh *m, double gradingFactor, double size[2], apf::Adjacent edgAdjVert, apf::Adjacent vertAdjEdg, std::queue< apf::MeshEntity * > &markedEdges, apf::MeshTag *isMarked, int fieldType, int vecPos, int idxFlag)
Definition: SizeField.cpp:1512
Smoother::apply
virtual void apply()
Definition: SizeField.cpp:1105
edgeWalkerInfo::L_local
double L_local
Definition: SizeField.cpp:184
Smoother::sum
apf::DynamicVector sum
Definition: SizeField.cpp:1127
Smoother::emptyTag
apf::MeshTag * emptyTag
Definition: SizeField.cpp:1129
checkForPropagation
int checkForPropagation(apf::Mesh *m, edgeWalkerInfo inputObject)
Definition: SizeField.cpp:191
BFS_propagation
int BFS_propagation(apf::Mesh *m, std::queue< edgeWalkerInfo > &markedVertices)
Definition: SizeField.cpp:224
MeshAdaptPUMIDrvr::numAdaptSteps
int numAdaptSteps
Definition: MeshAdaptPUMI.h:99
MeshAdaptPUMIDrvr::calculateSizeField
int calculateSizeField(double L_band)
Definition: SizeField.cpp:84
edgeWalkerInfo::vertex
apf::MeshEntity * vertex
Definition: SizeField.cpp:181
MeshAdaptPUMIDrvr::nAdapt
int nAdapt
Definition: MeshAdaptPUMI.h:93
SortingStruct::operator<
bool operator<(const SortingStruct &other) const
Definition: SizeField.cpp:927
MeshAdaptPUMIDrvr::logging_config
std::string logging_config
Definition: MeshAdaptPUMI.h:116
MeshAdaptPUMIDrvr::testIsotropicSizeField
int testIsotropicSizeField()
Definition: SizeField.cpp:1497
max
#define max(a, b)
Definition: analyticalSolutions.h:14
MeshAdaptPUMIDrvr::adapt_type_config
std::string adapt_type_config
Definition: MeshAdaptPUMI.h:115
edgeWalkerInfo::actualPosition
apf::Vector3 actualPosition
Definition: SizeField.cpp:182
Smoother::field
apf::Field * field
Definition: SizeField.cpp:1124
MeshAdaptPUMIDrvr::gradeMesh
int gradeMesh(double gradationFactor)
Definition: SizeField.cpp:1660
Smoother::nApplied
int nApplied
Definition: SizeField.cpp:1130
edgeWalkerInfo::initialRank
int initialRank
Definition: SizeField.cpp:188