HL_SignedDistTraverse.cc 7.5 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
//
// Software License for AMDiS
//
// Copyright (c) 2010 Dresden University of Technology 
// All rights reserved.
// Authors: Simon Vey, Thomas Witkowski et al.
//
// This file is part of AMDiS
//
// See also license.opensource.txt in the distribution.


13
14
15
#include "HL_SignedDistTraverse.h"
#include "VelocityExtFromVelocityField.h"

16
void HL_SignedDistTraverse::initializeBoundary() 
17
18
{
  FUNCNAME("HL_SignedDistTraverse::initializeBoundary()");
Thomas Witkowski's avatar
Thomas Witkowski committed
19
 
20
21
22
23
24
25
26
27
28
  // ===== 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();
29
  if (static_cast<int>(locInd.size()) < nBasFcts)
30
    locInd.resize(nBasFcts);
31
32
33
34
35
36
37
38
39
  
  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);
40
  } else {
41
42
43
44
45
46
    elInfo = stack.traverseFirst(feSpace->getMesh(),
				 -1, 
				 Mesh::CALL_LEAF_EL | 
				 Mesh::FILL_BOUND |
				 Mesh::FILL_COORDS);
  }
47

48
49
50
  while(elInfo) {

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

103

104
void HL_SignedDistTraverse::HL_updateIteration()
105
{
Thomas Witkowski's avatar
blub  
Thomas Witkowski committed
106
107
  FUNCNAME("HL_SignedDistTraverse::HL_updateIteration()");

108
109
  // ===== Create DOF vector for the last iteration step. =====
  if (sDOld_DOF)
110
111
112
    delete sDOld_DOF;

  sDOld_DOF = new DOFVector<double>(feSpace, "sDOld_DOF");
113
114
115
  sDOld_DOF->copy(const_cast<DOFVector<double> &>(*sD_DOF));
  
  // ===== Gauss-Seidel or Jacobi iteration ? =====
116
117
118
119
  if (GaussSeidelFlag)
    update_DOF = sD_DOF;  
  else
    update_DOF = sDOld_DOF;  
120
121
122
123
124
125
  
  // ===== Iteration loop: proceed until tolerance is reached. =====
  TraverseStack stack;
  ElInfo *elInfo;
  tol_reached = false;
  int itCntr = 0;
126
  while (!tol_reached && itCntr != maxIt) {   
127
128
129
    ++itCntr;
    tol_reached = true;
    
130
131
132
133
134
    // ===== 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) {
135
136
137
138
139
140
141
142
143
144
      HL_elementUpdate(elInfo);
      elInfo = stack.traverseNext(elInfo);
    }
    
    // ===== Is tolerance reached ? =====
    tol_reached = checkTol();
    
    sDOld_DOF->copy(const_cast<DOFVector<double> &>(*sD_DOF));
  }
  
145

Thomas Witkowski's avatar
blub  
Thomas Witkowski committed
146
  MSG("Calculation of signed distance function via mesh traverse iteration:\n");
147
  if (GaussSeidelFlag)
Thomas Witkowski's avatar
blub  
Thomas Witkowski committed
148
    MSG("\tGauss-Seidel iteration\n");
149
  else
Thomas Witkowski's avatar
blub  
Thomas Witkowski committed
150
    MSG("\tJacobi iteration\n");
151
  
Thomas Witkowski's avatar
blub  
Thomas Witkowski committed
152
  MSG("\tnumber of iterations needed: %d\n", itCntr);
153
154
}

155

156
void HL_SignedDistTraverse::HL_elementUpdate(ElInfo *elInfo) 
157
158
{
  FUNCNAME("HL_SignedDistTraverse::HL_elementUpdate()");
159
   
160
161
  // ===== Get global indices of vertices of element. =====
  const int nBasFcts = feSpace->getBasisFcts()->getNumber();
162
  if (static_cast<int>(locInd.size()) < nBasFcts)
163
    locInd.resize(nBasFcts);
164
  
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
165
166
167
  feSpace->getBasisFcts()->getLocalIndices(const_cast<Element *>(elInfo->getElement()),
					   const_cast<DOFAdmin *>(feSpace->getAdmin()),
					   locInd);
168
169
  
  // ===== Hopf-Lax element update for each vertex of element. =====
170
  for (int i = 0; i <= dim; i++) {
171
172
173
174
    
    // ===== Calculate update for non-boundary vertex. =====
    if ((*bound_DOF)[locInd[i]] < 1.e-15) {
      //save permutation of vertexes for calculation of the velocity
175
176
177
178
      if (velExt != NULL)
	velExt->setPermutation(i, 1);
	
      double update = calcElementUpdate(elInfo, i, &(locInd[0]));
179
180
181
182
183
184
185
186
      // ---> 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;
187
188
189
190
191

	// If Distance is corrected, calculate new velocity.
	if(velExt != NULL)
	  velExt->calcVelocity(&(locInd[0]), i);
	
192
193
194
195
196
	// ---> for test purposes: count number of calculated updates
	++setUpdate_Cntr;
	// ---> end: for test purposes
      }
    }
197
  } 
198
199
}

200

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

247

248
bool HL_SignedDistTraverse::checkTol()
249
250
251
252
{
  DOFVector<double>::Iterator it_sD(sD_DOF, USED_DOFS);
  DOFVector<double>::Iterator it_sDOld(sDOld_DOF, USED_DOFS);
  
253
254
  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)
255
256
257
258
      return false;
  
  return true;
}