HL_SignedDistTraverse.cc 7.47 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
  DegreeOfFreedom *locInd = new DegreeOfFreedom[nBasFcts];
18
19
20
21
22
23
24
25
26
  
  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);
27
  } else {
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
    elInfo = stack.traverseFirst(feSpace->getMesh(),
				 -1, 
				 Mesh::CALL_LEAF_EL | 
				 Mesh::FILL_BOUND |
				 Mesh::FILL_COORDS);
  }
  while(elInfo) {

    // Set elInfo in case velocity extension from velocity field is used.
    if (velExt && velExtType.isSet(VEL_EXT_FROM_VEL_FIELD)) {
      ((VelocityExtFromVelocityField *)velExt)->setElInfo(elInfo);
    }
    
    // Get local indices of vertices.
    feSpace->getBasisFcts()->getLocalIndices(
			     const_cast<Element *>(elInfo->getElement()),
			     const_cast<DOFAdmin *>(feSpace->getAdmin()),
			     locInd); 
    
    // Get element status.
    elStatus = elLS->createElementLevelSet(elInfo);
    
    // Is element cut by the interface ?
    if (elStatus == ElementLevelSet::LEVEL_SET_BOUNDARY) {
      
      // Reset element distance vector.
      for (int i=0; i<=dim; ++i) {
	distVec[i] = inftyValue;
      }
      
      // Mark all vertices as boundary vertices.
      for (int i=0; i<=dim; ++i) {
	(*bound_DOF)[locInd[i]] = 1.0;
      }
      
      // Calculate distance for all vertices.
      if (bndElDist->calcDistOnBoundaryElement(elInfo, distVec) !=
	  ElementLevelSet::LEVEL_SET_BOUNDARY) {
	ERROR_EXIT("error in distance calculation !\n");
      }
      else {
	
	// If distance is smaller, correct to new distance.
	for (int i=0; i<=dim; ++i) {
	  
	  // --> for test purposes:
	  if (distVec[i] > 1000) {
	    cout << "\nElement: " << elInfo->getElement()->getIndex() << ", Knoten: " << i << " , keine Randwertinitialisierung !\n";
	  }
	  // --> end: for test purposes
	  
	  if ((*sD_DOF)[locInd[i]] > distVec[i]) {
	    (*sD_DOF)[locInd[i]] = distVec[i];
	    //If Distance is corrected, calculate new velocity.
	    if(velExt != NULL)
	      {
		velExt->calcVelocityBoundary(locInd, i);
	      }
	  }
	}
      }
    }  // end of: elStatus == ElementLevelSet::LEVEL_SET_BOUNDARY
    
    elInfo = stack.traverseNext(elInfo);
  }  // end of: mesh traverse
  
94
  delete [] locInd;
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
}

void 
HL_SignedDistTraverse::HL_updateIteration()
{
  // ===== Create DOF vector for the last iteration step. =====
  if (sDOld_DOF)
    DELETE sDOld_DOF;
  sDOld_DOF = NEW DOFVector<double>(feSpace, "sDOld_DOF");
  sDOld_DOF->copy(const_cast<DOFVector<double> &>(*sD_DOF));
  
  // ===== Gauss-Seidel or Jacobi iteration ? =====
  if (GaussSeidelFlag) {
    update_DOF = sD_DOF;
  }
  else {
    update_DOF = sDOld_DOF;
  }
  
  // ===== Iteration loop: proceed until tolerance is reached. =====
  TraverseStack stack;
  ElInfo *elInfo;
  tol_reached = false;
  int itCntr = 0;
  while (!tol_reached && itCntr != maxIt) {
    
    ++itCntr;
    tol_reached = true;
    
    // ===== 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) {
      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";
  if (GaussSeidelFlag) {
    cout << "\tGauss-Seidel iteration\n";
  }
  else {
    cout << "\tJacobi iteration\n";
  }
  cout << "\tnumber of iterations needed: " << itCntr << "\n\n";
  
  return;
}

void 
HL_SignedDistTraverse::HL_elementUpdate(ElInfo *elInfo) 
{
  FUNCNAME("HL_SignedDistTraverse::HL_elementUpdate()");
  
  double update;
  
  // ===== Get global indices of vertices of element. =====
  const int nBasFcts = feSpace->getBasisFcts()->getNumber();
163
  DegreeOfFreedom *locInd = new DegreeOfFreedom[nBasFcts];
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
  
  feSpace->getBasisFcts()->getLocalIndices(
			   const_cast<Element *>(elInfo->getElement()),
			   const_cast<DOFAdmin *>(feSpace->getAdmin()),
			   locInd); 
  
  // ===== Hopf-Lax element update for each vertex of element. =====
  for (int i=0; i<=dim; ++i) {
    
    // ===== Calculate update for non-boundary vertex. =====
    if ((*bound_DOF)[locInd[i]] < 1.e-15) {
      //save permutation of vertexes for calculation of the velocity
      if(velExt != NULL)
	{
	  velExt->setPermutation(i, 1);
	}
      update = calcElementUpdate(elInfo, i, locInd);
      // ---> 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;
	///If Distance is corrected, calculate new velocity.
	    if(velExt != NULL)
	      {
		velExt->calcVelocity(locInd, i);
	      }
	// ---> for test purposes: count number of calculated updates
	++setUpdate_Cntr;
	// ---> end: for test purposes
      }
    }
  }
  
201
  delete [] locInd;
202
203
}

204
205
206
double HL_SignedDistTraverse::calcElementUpdate(ElInfo *elInfo, 
						int nXh,
						const DegreeOfFreedom *locInd) 
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
{
  double update;
  FixVec<WorldVector<double> *, VERTEX> elVert(dim, NO_INIT);
  FixVec<double, VERTEX> uhVal(dim, NO_INIT);
  
  // ===== Get local indices of element vertices (except xh). =====
  int nYh, nZh, nWh;
  
  nYh = (nXh + 1) % (dim+1);
  nZh = (nXh + 2) % (dim+1);
  if (dim == 3) nWh = (nXh + 3) % (dim+1);
  
  // ===== 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) {
  case 2: elVert[0] = &(elInfo->getCoord(nYh));
    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;
  case 3: elVert[0] = &(elInfo->getCoord(nYh));
    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. =====
  update = elUpdate->calcElementUpdate(elVert, uhVal);
  
  return update;
}

254
bool HL_SignedDistTraverse::checkTol()
255
256
257
258
259
260
261
262
263
264
265
266
267
268
{
  DOFVector<double>::Iterator it_sD(sD_DOF, USED_DOFS);
  DOFVector<double>::Iterator it_sDOld(sDOld_DOF, USED_DOFS);
  
  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) {
      return false;
    }
  }
  
  return true;
}