HL_SignedDistTraverse.cc 7.22 KB
Newer Older
1
2
3
#include "HL_SignedDistTraverse.h"
#include "VelocityExtFromVelocityField.h"

4
void HL_SignedDistTraverse::initializeBoundary() 
5
6
7
8
9
10
11
12
13
14
15
16
{
  FUNCNAME("HL_SignedDistTraverse::initializeBoundary()");
  
  // ===== All non-boundary vertices are initialized with "infinity". =====
  sD_DOF->set(inftyValue);
  
  // ===== Traverse mesh and initialize boundary elements. =====
  TraverseStack stack;
  FixVec<double, VERTEX> distVec(dim, NO_INIT);
  int elStatus;
  
  const int nBasFcts = feSpace->getBasisFcts()->getNumber();
17
18
  if (locInd.size() < nBasFcts)
    locInd.resize(nBasFcts);
19
20
21
22
23
24
25
26
27
  
  ElInfo *elInfo;
  if (velExt && velExtType.isSet(VEL_EXT_FROM_VEL_FIELD)) {
    elInfo = stack.traverseFirst(feSpace->getMesh(),
				 -1, 
				 Mesh::CALL_LEAF_EL | 
				 Mesh::FILL_BOUND |
				 Mesh::FILL_COORDS |
				 Mesh::FILL_GRD_LAMBDA);
28
  } else {
29
30
31
32
33
34
    elInfo = stack.traverseFirst(feSpace->getMesh(),
				 -1, 
				 Mesh::CALL_LEAF_EL | 
				 Mesh::FILL_BOUND |
				 Mesh::FILL_COORDS);
  }
35

36
37
38
  while(elInfo) {

    // Set elInfo in case velocity extension from velocity field is used.
39
40
    if (velExt && velExtType.isSet(VEL_EXT_FROM_VEL_FIELD))
      ((VelocityExtFromVelocityField *)velExt)->setElInfo(elInfo);    
41
42
43
    
    // Get local indices of vertices.
    feSpace->getBasisFcts()->getLocalIndices(
44
45
46
			     const_cast<Element*>(elInfo->getElement()),
			     const_cast<DOFAdmin*>(feSpace->getAdmin()),
			     &(locInd[0])); 
47
48
49
50
51
52
53
54
    
    // Get element status.
    elStatus = elLS->createElementLevelSet(elInfo);
    
    // Is element cut by the interface ?
    if (elStatus == ElementLevelSet::LEVEL_SET_BOUNDARY) {
      
      // Reset element distance vector.
55
56
      for (int i=0; i <= dim; i++)
	distVec[i] = inftyValue;      
57
58
      
      // Mark all vertices as boundary vertices.
59
60
      for (int i = 0; i <= dim; i++)
	(*bound_DOF)[locInd[i]] = 1.0;      
61
62
63
64
65
      
      // Calculate distance for all vertices.
      if (bndElDist->calcDistOnBoundaryElement(elInfo, distVec) !=
	  ElementLevelSet::LEVEL_SET_BOUNDARY) {
	ERROR_EXIT("error in distance calculation !\n");
66
      } else {
67
68
	
	// If distance is smaller, correct to new distance.
69
	for (int i = 0; i <= dim; i++) {
70
71
	  
	  // --> for test purposes:
72
73
74
75
	  if (distVec[i] > 1000)
	    cout << "\nElement: " << elInfo->getElement()->getIndex() 
		 << ", Knoten: " << i << " , keine Randwertinitialisierung !\n";
	  
76
77
78
79
80
	  // --> end: for test purposes
	  
	  if ((*sD_DOF)[locInd[i]] > distVec[i]) {
	    (*sD_DOF)[locInd[i]] = distVec[i];
	    //If Distance is corrected, calculate new velocity.
81
82
	    if (velExt != NULL)	      
	      velExt->calcVelocityBoundary(&(locInd[0]), i);	      
83
84
85
86
87
88
	  }
	}
      }
    }  // end of: elStatus == ElementLevelSet::LEVEL_SET_BOUNDARY
    
    elInfo = stack.traverseNext(elInfo);
89
  }  // end of: mesh traverse 
90
91
}

92
void HL_SignedDistTraverse::HL_updateIteration()
93
94
95
{
  // ===== Create DOF vector for the last iteration step. =====
  if (sDOld_DOF)
96
97
98
    delete sDOld_DOF;

  sDOld_DOF = new DOFVector<double>(feSpace, "sDOld_DOF");
99
100
101
  sDOld_DOF->copy(const_cast<DOFVector<double> &>(*sD_DOF));
  
  // ===== Gauss-Seidel or Jacobi iteration ? =====
102
103
104
105
  if (GaussSeidelFlag)
    update_DOF = sD_DOF;  
  else
    update_DOF = sDOld_DOF;  
106
107
108
109
110
111
  
  // ===== Iteration loop: proceed until tolerance is reached. =====
  TraverseStack stack;
  ElInfo *elInfo;
  tol_reached = false;
  int itCntr = 0;
112
  while (!tol_reached && itCntr != maxIt) {   
113
114
115
    ++itCntr;
    tol_reached = true;
    
116
117
118
119
120
    // ===== Traverse mesh: perform Hopf-Lax element update on each element. =====
    elInfo = 
      stack.traverseFirst(feSpace->getMesh(), -1, 
			  Mesh::CALL_LEAF_EL | Mesh::FILL_BOUND | Mesh::FILL_COORDS);
    while (elInfo) {
121
122
123
124
125
126
127
128
129
130
131
      HL_elementUpdate(elInfo);
      elInfo = stack.traverseNext(elInfo);
    }
    
    // ===== Is tolerance reached ? =====
    tol_reached = checkTol();
    
    sDOld_DOF->copy(const_cast<DOFVector<double> &>(*sD_DOF));
  }
  
  cout << "\n\nCalculation of signed distance function via mesh traverse iteration:\n";
132
133
134
135

  if (GaussSeidelFlag)
    cout << "\tGauss-Seidel iteration\n";  
  else
136
137
    cout << "\tJacobi iteration\n";
  
138
  cout << "\tnumber of iterations needed: " << itCntr << "\n\n";
139
140
}

141
void HL_SignedDistTraverse::HL_elementUpdate(ElInfo *elInfo) 
142
143
{
  FUNCNAME("HL_SignedDistTraverse::HL_elementUpdate()");
144
   
145
146
  // ===== Get global indices of vertices of element. =====
  const int nBasFcts = feSpace->getBasisFcts()->getNumber();
147
148
  if (locInd.size() < nBasFcts)
    locInd.resize(nBasFcts);
149
150
151
152
  
  feSpace->getBasisFcts()->getLocalIndices(
			   const_cast<Element *>(elInfo->getElement()),
			   const_cast<DOFAdmin *>(feSpace->getAdmin()),
153
			   &(locInd[0])); 
154
155
  
  // ===== Hopf-Lax element update for each vertex of element. =====
156
  for (int i = 0; i <= dim; i++) {
157
158
159
160
    
    // ===== Calculate update for non-boundary vertex. =====
    if ((*bound_DOF)[locInd[i]] < 1.e-15) {
      //save permutation of vertexes for calculation of the velocity
161
162
163
164
      if (velExt != NULL)
	velExt->setPermutation(i, 1);
	
      double update = calcElementUpdate(elInfo, i, &(locInd[0]));
165
166
167
168
169
170
171
172
      // ---> for test purposes: count number of calculated updates
      ++calcUpdate_Cntr;
      // ---> end: for test purposes
      
      // Calculates minimum of all element updates on elements
      // containing vertex i.
      if (update < (*sD_DOF)[locInd[i]]) {
	(*sD_DOF)[locInd[i]] = update;
173
174
175
176
177

	// If Distance is corrected, calculate new velocity.
	if(velExt != NULL)
	  velExt->calcVelocity(&(locInd[0]), i);
	
178
179
180
181
182
	// ---> for test purposes: count number of calculated updates
	++setUpdate_Cntr;
	// ---> end: for test purposes
      }
    }
183
  } 
184
185
}

186
187
188
double HL_SignedDistTraverse::calcElementUpdate(ElInfo *elInfo, 
						int nXh,
						const DegreeOfFreedom *locInd) 
189
190
{
  // ===== Get local indices of element vertices (except xh). =====
191
192
193
194
195
  int nWh = 0;  
  int nYh = (nXh + 1) % (dim+1);
  int nZh = (nXh + 2) % (dim+1);
  if (dim == 3) 
    nWh = (nXh + 3) % (dim+1);
196
197
198
199
200
201
  
  // ===== Get world coordinates of vertices of element and their values
  //       of uh. 
  //       The coordinates of the vertex the update is calculated for
  //       are stored at the end of the vector elVert. =====
  switch (dim) {
202
203
  case 2: 
    elVert[0] = &(elInfo->getCoord(nYh));
204
205
206
207
208
209
210
211
    elVert[1] = &(elInfo->getCoord(nZh));
    elVert[2] = &(elInfo->getCoord(nXh));
    
    uhVal[0] = (*update_DOF)[locInd[nYh]];
    uhVal[1] = (*update_DOF)[locInd[nZh]];
    uhVal[2] = (*update_DOF)[locInd[nXh]];
    
    break;
212
213
  case 3: 
    elVert[0] = &(elInfo->getCoord(nYh));
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
    elVert[1] = &(elInfo->getCoord(nZh));
    elVert[2] = &(elInfo->getCoord(nWh));
    elVert[3] = &(elInfo->getCoord(nXh));

    uhVal[0] = (*update_DOF)[locInd[nYh]];
    uhVal[1] = (*update_DOF)[locInd[nZh]];
    uhVal[2] = (*update_DOF)[locInd[nWh]];
    uhVal[3] = (*update_DOF)[locInd[nXh]];
    
    break;
  default: ERROR_EXIT("illegal dimension !\n");
    break;
  }
  
  // ===== Calculate Hopf-Lax element update for vertex. =====
229
  return elUpdate->calcElementUpdate(elVert, uhVal);
230
231
}

232
bool HL_SignedDistTraverse::checkTol()
233
234
235
236
{
  DOFVector<double>::Iterator it_sD(sD_DOF, USED_DOFS);
  DOFVector<double>::Iterator it_sDOld(sDOld_DOF, USED_DOFS);
  
237
238
  for (it_sD.reset(), it_sDOld.reset(); !it_sD.end(); ++it_sD, ++it_sDOld)
    if ((*it_sDOld) - (*it_sD) > tol || (*it_sDOld) - (*it_sD) < 0)
239
240
241
242
      return false;
  
  return true;
}