Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
turbulent_viscosity_fvcr_impl.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2013: G-CSC, Goethe University Frankfurt
3 * Author: Christian Wehner
4 *
5 * This file is part of UG4.
6 *
7 * UG4 is free software: you can redistribute it and/or modify it under the
8 * terms of the GNU Lesser General Public License version 3 (as published by the
9 * Free Software Foundation) with the following additional attribution
10 * requirements (according to LGPL/GPL v3 §7):
11 *
12 * (1) The following notice must be displayed in the Appropriate Legal Notices
13 * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 *
15 * (2) The following notice must be displayed at a prominent place in the
16 * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 *
18 * (3) The following bibliography is recommended for citation and must be
19 * preserved in all covered files:
20 * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 * parallel geometric multigrid solver on hierarchically distributed grids.
22 * Computing and visualization in science 16, 4 (2013), 151-164"
23 * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 * flexible software system for simulating pde based models on high performance
25 * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 *
27 * This program is distributed in the hope that it will be useful,
28 * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 * GNU Lesser General Public License for more details.
31 */
32
33#ifndef __H__UG__NAVIER_STOKES_TURBULENT_VISCOSITY_DATA_IMPL_H_
34#define __H__UG__NAVIER_STOKES_TURBULENT_VISCOSITY_DATA_IMPL_H_
35
36namespace ug{
37namespace NavierStokes{
38
39template <typename TData, int dim, typename TImpl,typename TGridFunction>
41 for(int si = 0; si < approximationSpace.domain()->subset_handler()->num_subsets(); ++si){
42 // transfer to lower levels, averaging over child edges (2d) / child faces (3d)
43 for (size_t lev=approximationSpace.num_levels()-2;(int)lev>=0;lev--){
44 const DoFDistribution& lDD = *approximationSpace.dof_distribution(GridLevel(lev, GridLevel::LEVEL));
45 const MultiGrid& grid = *lDD.multi_grid();
46 typedef typename DoFDistribution::traits<side_type>::const_iterator coarseLevelSideIter;
47 coarseLevelSideIter clsIter, clsIterEnd;
48 clsIter = lDD.template begin<side_type>(si);
49 clsIterEnd = lDD.template end<side_type>(si);
50 for (;clsIter != clsIterEnd;clsIter++){
51 side_type* elem = *clsIter;
52 size_t numChildren = grid.num_children<side_type>(elem);
53 number avgValue=0;
54 for (size_t i=0;i<numChildren;i++){
55 avgValue += aaData[grid.get_child<side_type>(elem, i)];
56 }
57 avgValue/=numChildren;
58 aaData[elem] = avgValue;
59 }
60 if (lev==0) break;
61 }
62 }
63}
64
65template <typename TData, int dim, typename TImpl,typename TGridFunction>
67 // get domain
68 domain_type& domain = *u->domain().get();
69 // create Multiindex
70 std::vector<DoFIndex> multInd;
71 for(int si = 0; si < domain.subset_handler()->num_subsets(); ++si){
72 ElemIterator iter = u->template begin<side_type>(si);
73 ElemIterator iterEnd = u->template end<side_type>(si);
74 for( ;iter !=iterEnd; ++iter)
75 {
76 side_type* side = iter;
77 for (int d=0;d<dim;d++){
78 u->dof_indices(side, d, multInd);
79 aaU[side][d]=DoFRef(*u,multInd[0]);
80 }
81 }
82 }
83}
84
85void assignValue(number& v,size_t ind,number value){
86 v=value;
87}
88
89void assignValue(MathVector<2>& v,size_t ind,number value){
90 v[ind]=value;
91}
92
93void assignValue(MathVector<3>& v,size_t ind,number value){
94 v[ind]=value;
95}
96
97void assignValue(MathSymmetricMatrix<2>& v,size_t ind,number value){
98 v[ind]=value;
99}
100
101void assignValue(MathSymmetricMatrix<3>& v,size_t ind,number value){
102 v[ind]=value;
103}
104
105// go over all elements, interpolate data to barycenter, average by multiplying with corresponding element volume and deviding by complete adjacent element volume
106// parameters: filtered values, filter volume (computed in function), original values
107template <typename TData, int dim, typename TImpl,typename TGridFunction>
108template <typename VType>
111 bool useGridFunction = true;
112 if (u==SPNULL) useGridFunction = false;
113
114 std::vector<DoFIndex> multInd;
115
116 // get domain of grid function
117 domain_type& domain = *m_uInfo->domain().get();
119
120 // set attachment values to zero
121 SetAttachmentValues(aaUHat , m_uInfo->template begin<side_type>(), m_uInfo->template end<side_type>(), 0);
122 SetAttachmentValues(aaVol , m_uInfo->template begin<side_type>(), m_uInfo->template end<side_type>(), 0);
123
124 // coord and vertex array
127
128 // get position accessor
129 typedef typename domain_type::position_accessor_type position_accessor_type;
130 const position_accessor_type& posAcc = domain.position_accessor();
131
132 // assemble deformation tensor fluxes
133 for(int si = 0; si < domain.subset_handler()->num_subsets(); ++si)
134 {
135 // get iterators
136 ElemIterator iter = m_uInfo->template begin<elem_type>(si);
137 ElemIterator iterEnd = m_uInfo->template end<elem_type>(si);
138
139 // loop elements of dimension
140 for( ;iter !=iterEnd; ++iter)
141 {
142 // get Elem
143 elem_type* elem = *iter;
144 // get vertices and extract corner coordinates
145 const size_t numVertices = elem->num_vertices();
146 MathVector<dim> bary,localBary;
147 bary = 0;
148
149 for(size_t i = 0; i < numVertices; ++i){
150 vVrt[i] = elem->vertex(i);
151 coCoord[i] = posAcc[vVrt[i]];
152 bary += coCoord[i];
153 };
154 bary /= numVertices;
155
156 // reference object id
157 ReferenceObjectID roid = elem->reference_object_id();
158
159 // get trial space
160 const LocalShapeFunctionSet<dim>& rTrialSpace =
161 LocalFiniteElementProvider::get<dim>(roid, LFEID(LFEID::CROUZEIX_RAVIART, dim, 1));
162
163 // get Reference Mapping
164 DimReferenceMapping<dim, dim>& map = ReferenceMappingProvider::get<dim, dim>(roid, coCoord);
165
166 map.global_to_local(localBary,bary);
167
168 typedef typename grid_type::template traits<side_type>::secure_container side_secure_container;
169 side_secure_container sides;
170
171 UG_ASSERT(dynamic_cast<elem_type*>(elem) != NULL, "Only elements of type elem_type are currently supported");
172
173 domain.grid()->associated_elements_sorted(sides, static_cast<elem_type*>(elem) );
174 if (m_bAdaptive){
175 std::vector<DoFIndex> ind;
176 m_uInfo->dof_indices(elem,0,ind,true,true);
177 get_constrained_sides_cr<side_type,side_secure_container,TGridFunction>(sides,*m_uInfo,ind);
178 }
179
180 // memory for shapes
181 std::vector<number> vShape;
182
183 // evaluate finite volume geometry
184 if (!m_bAdaptive)
185 geo.update(elem, &(coCoord[0]), domain.subset_handler().get());
186 else
187 geo.update_hanging(elem, &(coCoord[0]), domain.subset_handler().get());
188
189 rTrialSpace.shapes(vShape, localBary);
190
191 size_t nofsides = geo.num_scv();
192
193 VType value;
194 value = 0;
195 number elementVolume = 0;
196 for (size_t i=0;i<nofsides;i++){
197 const typename DimCRFVGeometry<dim>::SCV& scv = geo.scv(i);
198 size_t s = scv.node_id();
199 VType localValue;
200 if (useGridFunction)
201 for (int d=0;d<dim;d++){
202 u->dof_indices(sides[s], d, multInd);
203 assignValue(localValue,d,DoFRef(*u,multInd[0]));
204 }
205 else
206 localValue = (*aaU)[sides[s]];
207 localValue *= vShape[s];
208 value += localValue;
209 elementVolume += scv.volume();
210 }
211 value *= elementVolume;
212 for (size_t s=0;s<nofsides;s++){
213 aaUHat[sides[s]] += value;
214 aaVol[sides[s]] += elementVolume;
215 }
216 }
217 }
218 PeriodicBoundaryManager* pbm = (domain.grid())->periodic_boundary_manager();
219 // average
220 for(int si = 0; si < domain.subset_handler()->num_subsets(); ++si)
222 SideIterator sideIter = m_uInfo->template begin<side_type>(si);
223 SideIterator sideIterEnd = m_uInfo->template end<side_type>(si);
224 // no filtering in boundary nodes
225 if ((m_turbZeroSg.size()!=0) && (m_turbZeroSg.contains(si)==true)){
226 for( ;sideIter !=sideIterEnd; sideIter++)
227 {
228 side_type* side = *sideIter;
229 if (pbm && pbm->is_slave(side)) continue;
230 if (useGridFunction){
231 for (int d=0;d<dim;d++){
232 u->dof_indices(side, d, multInd);
233 assignValue(aaUHat[side],d,DoFRef(*u,multInd[0]));
234 }
236 else
237 aaUHat[side] = (*aaU)[side];
238 }
239 } else {
240 for( ;sideIter !=sideIterEnd; sideIter++)
241 {
242 side_type* side = *sideIter;
243 if (pbm && pbm->is_slave(side)) continue;
244 // constrained edges
245 if (aaVol[side]==0) continue;
246 aaUHat[side]/=(number)aaVol[side];
247 }
248 }
249 }
250 if (m_bAdaptive){
251 if (dim==2) constrainingSideAveraging<TGridFunction,side_type,ConstrainingEdge,VType>(aaUHat,m_uInfo);
252 else {
253 constrainingSideAveraging<TGridFunction,side_type,ConstrainingTriangle,VType>(aaUHat,m_uInfo);
254 constrainingSideAveraging<TGridFunction,side_type,ConstrainingQuadrilateral,VType>(aaUHat,m_uInfo);
255 }
256 }
257}
258
259// go over all elements, interpolate data to scv barycenter, average by multiplying with corresponding scv volume and deviding by volume of complete control volume
260template <typename TData, int dim, typename TImpl,typename TGridFunction>
261template <typename VType>
264 bool useGridFunction = true;
265 if (u==NULL) useGridFunction = false;
266
267 std::vector<DoFIndex> multInd;
268 // get domain of grid function
269 domain_type& domain = *m_uInfo->domain().get();
271
272 // set attachment values to zero
273 SetAttachmentValues(aaUHat , m_uInfo->template begin<side_type>(), m_uInfo->template end<side_type>(), 0);
274 SetAttachmentValues(aaVol , m_uInfo->template begin<side_type>(), m_uInfo->template end<side_type>(), 0);
275
276 // coord and vertex array
279
280 // get position accessor
281 typedef typename domain_type::position_accessor_type position_accessor_type;
282 const position_accessor_type& posAcc = domain.position_accessor();
283
284 // assemble deformation tensor fluxes
285 for(int si = 0; si < domain.subset_handler()->num_subsets(); ++si)
286 {
287 // get iterators
288 ElemIterator iter = m_uInfo->template begin<elem_type>(si);
289 ElemIterator iterEnd = m_uInfo->template end<elem_type>(si);
290
291 // loop elements of dimension
292 for( ;iter !=iterEnd; ++iter)
293 {
294 // get Elem
295 elem_type* elem = *iter;
296 // get vertices and extract corner coordinates
297 const size_t numVertices = elem->num_vertices();
298 MathVector<dim> bary,localBary;
299 bary = 0;
300
301 for(size_t i = 0; i < numVertices; ++i){
302 vVrt[i] = elem->vertex(i);
303 coCoord[i] = posAcc[vVrt[i]];
304 bary += coCoord[i];
305 };
306 bary /= numVertices;
307
308 // reference object id
309 ReferenceObjectID roid = elem->reference_object_id();
310
311 // get trial space
312 const LocalShapeFunctionSet<dim>& rTrialSpace =
313 LocalFiniteElementProvider::get<dim>(roid, LFEID(LFEID::CROUZEIX_RAVIART, dim, 1));
314
315 // get Reference Mapping
316 DimReferenceMapping<dim, dim>& map = ReferenceMappingProvider::get<dim, dim>(roid, coCoord);
317
318 map.global_to_local(localBary,bary);
319
320 typedef typename grid_type::template traits<side_type>::secure_container side_secure_container;
321 side_secure_container sides;
322
323 UG_ASSERT(dynamic_cast<elem_type*>(elem) != NULL, "Only elements of type elem_type are currently supported");
324
325 domain.grid()->associated_elements_sorted(sides, static_cast<elem_type*>(elem) );
326
327 // evaluate finite volume geometry
328 if (!m_bAdaptive){
329 geo.update(elem, &(coCoord[0]), domain.subset_handler().get());
330 } else {
331 geo.update_hanging(elem, &(coCoord[0]), domain.subset_handler().get());
332 if (geo.num_constrained_dofs()>0){
333 m_uInfo->dof_indices(elem,0,multInd,true,true);
334 get_constrained_sides_cr<side_type,side_secure_container,TGridFunction>(sides,*m_uInfo,multInd);
335 }
336 }
337
338 size_t nNaturalSides = geo.num_sh();
339 size_t nSides = geo.num_scv();
340
341 MathVector<dim> scvBary,scvLocalBary;
342 for (size_t i=0;i<nSides;i++){
343 const typename DimCRFVGeometry<dim>::SCV& scv = geo.scv(i);
344 size_t s = geo.node_id();
345 scvBary = 0;
346 // compute barycenter of scv (average of side corner nodes + element bary)
347 for (size_t j=0;j<scv.num_corners();j++){
348 scvBary += scv.global_corner(j);
349 }
350 scvBary/=(number)scv.num_corners();
351 map.global_to_local(scvLocalBary,scvBary);
352 // memory for shapes
353 std::vector<number> vShape;
354 rTrialSpace.shapes(vShape, scvLocalBary);
355 VType localValue = 0;
356 for (size_t j=0;j<nNaturalSides;j++){
357 if (useGridFunction){
358 VType uValue;
359 for (int d=0;d<dim;d++){
360 u->dof_indices(sides[j], d, multInd);
361 assignValue(uValue,d,DoFRef(*u,multInd[0]));
362 }
363 localValue += vShape[j]*uValue;
364 } else
365 localValue += vShape[j]*aaU[sides[j]];
366 }
367 localValue *= scv.volume();
368 aaVol[sides[s]] += scv.volume();
369 aaUHat[sides[s]] += localValue;
370 }
371 }
372 }
373 PeriodicBoundaryManager* pbm = (domain.grid())->periodic_boundary_manager();
374 // average
375 for(int si = 0; si < domain.subset_handler()->num_subsets(); ++si)
376 {
377 SideIterator sideIter = m_uInfo->template begin<side_type>(si);
378 SideIterator sideIterEnd = m_uInfo->template end<side_type>(si);
379 // no filtering in boundary nodes
380 if ((m_turbZeroSg.size()!=0) && (m_turbZeroSg.contains(si)==true)){
381 for( ;sideIter !=sideIterEnd; sideIter++)
382 {
383 side_type* side = *sideIter;
384 if (pbm && pbm->is_slave(side)) continue;
385 aaUHat[side] = aaU[side];
386 }
387 } else {
388 for( ;sideIter !=sideIterEnd; sideIter++)
389 {
390 side_type* side = *sideIter;
391 if (pbm && pbm->is_slave(side)) continue;
392 aaUHat[side]/=(number)aaVol[side];
393 }
394 }
395 }
396 if (m_bAdaptive){
397 if (dim==2) constrainingSideAveraging<TGridFunction,side_type,ConstrainingEdge,VType>(aaUHat,m_uInfo);
398 else {
399 constrainingSideAveraging<TGridFunction,side_type,ConstrainingTriangle,VType>(aaUHat,m_uInfo);
400 constrainingSideAveraging<TGridFunction,side_type,ConstrainingQuadrilateral,VType>(aaUHat,m_uInfo);
401 }
402 }
403}
404
405template <typename TData, int dim, typename TImpl,typename TGridFunction>
407 bool useGridFunction = true;
408 if (u==SPNULL) useGridFunction = false;
409
410 // get domain
411 domain_type& domain = *m_uInfo->domain().get();
412 // get periodic boundary manager
413 PeriodicBoundaryManager* pbm = (domain.grid())->periodic_boundary_manager();
414
415 // create Multiindex
416 std::vector<DoFIndex> multInd;
417
419
420 // add boundary subsets to enforce boundary subset computations in geo.update()
421 for(size_t i = 0; i < m_turbZeroSg.size(); ++i){
422 geo.add_boundary_subset(m_turbZeroSg[i]);
423 }
424
425 // set attachment values to zero
426 SetAttachmentValues(aaDefTensor , m_uInfo->template begin<side_type>(), m_uInfo->template end<side_type>(), 0);
427 SetAttachmentValues(aaVol , m_uInfo->template begin<side_type>(), m_uInfo->template end<side_type>(), 0);
428
429 // coord and vertex array
432
433 // get position accessor
434 typedef typename domain_type::position_accessor_type position_accessor_type;
435 const position_accessor_type& posAcc = domain.position_accessor();
436
437 for(int si = 0; si < domain.subset_handler()->num_subsets(); ++si){
438 // get iterators
439 ElemIterator iter = m_uInfo->template begin<elem_type>(si);
440 ElemIterator iterEnd = m_uInfo->template end<elem_type>(si);
441
442 // loop elements of dimension
443 for( ;iter !=iterEnd; ++iter)
444 {
445 // get Elem
446 elem_type* elem = *iter;
447 // get vertices and extract corner coordinates
448 const size_t numVertices = elem->num_vertices();
449 for(size_t i = 0; i < numVertices; ++i){
450 vVrt[i] = elem->vertex(i);
451 coCoord[i] = posAcc[vVrt[i]];
452 };
453
454 static const size_t MaxNumSidesOfElem = 10;
455
456 typedef MathVector<dim> MVD;
457 std::vector<MVD> uValue(MaxNumSidesOfElem);
458 MVD ipVelocity;
459
460 typedef typename grid_type::template traits<side_type>::secure_container side_secure_container;
461 side_secure_container sides;
462
463 UG_ASSERT(dynamic_cast<elem_type*>(elem) != NULL, "Only elements of type elem_type are currently supported");
464
465 domain.grid()->associated_elements_sorted(sides, static_cast<elem_type*>(elem) );
466
467 // evaluate finite volume geometry
468 if (!m_bAdaptive){
469 geo.update(elem, &(coCoord[0]), domain.subset_handler().get());
470 } else {
471 geo.update_hanging(elem, &(coCoord[0]), domain.subset_handler().get());
472 if (geo.num_constrained_dofs()>0){
473 m_uInfo->dof_indices(elem,0,multInd,true,true);
474 get_constrained_sides_cr<side_type,side_secure_container,TGridFunction>(sides,*m_uInfo,multInd);
475 }
476 }
477
478 size_t nNaturalSides = geo.num_sh();
479 size_t nSides = geo.num_scv();
480 size_t nip = geo.num_scvf();
481
482 if (useGridFunction){
483 for (size_t s=0;s < nNaturalSides;s++){
484 for (int d=0;d<dim;d++){
485 m_uInfo->dof_indices(sides[s], d, multInd);
486 uValue[s][d]=DoFRef(*u,multInd[0]);
487 }
488 }
489 } else
490 {
491 for (size_t s=0;s < nNaturalSides;s++)
492 uValue[s]=(*aaU)[sides[s]];
493 }
494
495 for (size_t i=0;i < nSides;i++)
496 {
497 const typename DimCRFVGeometry<dim>::SCV& scv = geo.scv(i);
498 size_t s = scv.node_id();
499 aaVol[sides[s]] += scv.volume();
500 }
501
502 for (size_t ip=0;ip<nip;ip++){
503 // get current SCVF
504 ipVelocity = 0;
505 const typename DimCRFVGeometry<dim>::SCVF& scvf = geo.scvf(ip);
506 for (size_t s=0;s < nNaturalSides;s++){
507 for (int d=0;d<dim;d++){
508 ipVelocity[d] += scvf.shape(s)*uValue[s][d];
509 };
510 };
511 dimMat ipDefTensorFlux;
512 ipDefTensorFlux = 0;
513 size_t ind = 0;
514 for (int i=0;i<dim;i++){
515 for (int j=i;j<dim;j++){
516 ipDefTensorFlux[ind]= 0.5 * (ipVelocity[i] * scvf.normal()[j] + ipVelocity[j] * scvf.normal()[i]);
517 ind++;
518 }
519 }
520 aaDefTensor[sides[scvf.from()]]+=ipDefTensorFlux;
521 aaDefTensor[sides[scvf.to()]]-=ipDefTensorFlux;
522 }
523 for(size_t sgi = 0; sgi < m_turbZeroSg.size(); ++sgi){
524 const size_t sgsi=m_turbZeroSg[sgi];
525 if (geo.num_bf(sgsi) == 0) continue;
526 for(size_t bfi = 0; bfi < geo.num_bf(sgsi); ++bfi){
527 const typename DimCRFVGeometry<dim>::BF& bf = geo.bf(sgsi, bfi);
528 const size_t sideID = bf.node_id();
529 size_t ind = 0;
530 for (int i=0;i<dim;i++){
531 for (int j=i;j<dim;j++){
532 // bf ip and u position are identical for CR-FV-Geometry
533 //for debug UG_LOG("[" << i << "," << j << "]" << 0.5 * (uValue[sideID][i] * bf.normal()[j] + uValue[sideID][j] * bf.normal()[i]) << "\n");
534 aaDefTensor[sides[sideID]][ind] -= 0.5 * (uValue[sideID][i] * bf.normal()[j] + uValue[sideID][j] * bf.normal()[i]);
535 ind++;
536 }
537 }
538 }
539 }
540 }
541 }
542 // average
543 for(int si = 0; si < domain.subset_handler()->num_subsets(); ++si){
544 SideIterator sideIter = m_uInfo->template begin<side_type>(si);
545 SideIterator sideIterEnd = m_uInfo->template end<side_type>(si);
546 //for debug UG_LOG("|||||||||||||||||||||||||||||||||||\n");
547 //for debug UG_LOG("si=" << si << "\n");
548 for( ;sideIter !=sideIterEnd; sideIter++)
549 {
550 side_type* side = *sideIter;
551 if (pbm && pbm->is_slave(side)) continue;
552 aaDefTensor[side]/=(number)aaVol[side];
553 //for debug UG_LOG("$$$$$$$$$$$$$$$$$$$$$$$$$$$$\n");
554 // MathVector<dim> posCo;
555 // for (int d=0;d<dim;d++)
556 // posCo[d] = 0.5*posAcc[side->vertex(0)][d] + 0.5*posAcc[side->vertex(1)][d];
557 //for debug UG_LOG(" c=" << posCo << "\n");
558 //for debug for (int d1=0;d1<dim;d1++)
559 //for debug for (int d2=0;d2<dim;d2++)
560 //for debug UG_LOG(" tensor(" << d1 << "," << d2 << ")=" << aaDefTensor[side][d1][d2] << "\n");
561 //for debug UG_LOG(" norm=" << FNorm(aaDefTensor[side]) << "\n");
562 }
563 }
564 // handle constrained sides
565 if (m_bAdaptive){
566 if (dim==2) constrainingSideAveraging<TGridFunction,side_type,ConstrainingEdge,dimMat>(aaDefTensor,m_uInfo);
567 else {
568 constrainingSideAveraging<TGridFunction,side_type,ConstrainingTriangle,dimMat>(aaDefTensor,m_uInfo);
569 constrainingSideAveraging<TGridFunction,side_type,ConstrainingQuadrilateral,dimMat>(aaDefTensor,m_uInfo);
570 }
571 }
572}
573
574template <typename TData, int dim, typename TImpl,typename TGridFunction>
576 // get domain of grid function
577 domain_type& domain = *m_uInfo->domain().get();
578 // get periodic boundary manager
579 PeriodicBoundaryManager* pbm = (domain.grid())->periodic_boundary_manager();
580 // average
581 for(int si = 0; si < domain.subset_handler()->num_subsets(); ++si){
582 SideIterator sideIter = m_uInfo->template begin<side_type>(si);
583 SideIterator sideIterEnd = m_uInfo->template end<side_type>(si);
584 for( ;sideIter !=sideIterEnd; sideIter++)
585 {
586 side_type* side = *sideIter;
587 if (pbm && pbm->is_slave(side)) continue;
588 //UG_LOG(aaTensor[side] << "\n");
589 //UG_LOG(FNorm(aaTensor[side]) << "\n");
590 //UG_LOG("--------------------\n");
591 //for debug UG_LOG("&&&&&&&& norm = " << FNorm(aaTensor[side]) << "\n");
592 aaTensor[side].scale_by_fnorm();
593 }
594 }
595}
596
597template <typename TData, int dim, typename TImpl,typename TGridFunction>
600 bool useGridFunction = true;
601 if (u==SPNULL) useGridFunction = false;
602
603 // get domain of grid function
604 domain_type& domain = *m_uInfo->domain().get();
605 // get periodic boundary manager
606 PeriodicBoundaryManager* pbm = (domain.grid())->periodic_boundary_manager();
607 // create Multiindex
608 std::vector<DoFIndex> multInd;
609
610 for(int si = 0; si < domain.subset_handler()->num_subsets(); ++si){
611 SideIterator sideIter = m_uInfo->template begin<side_type>(si);
612 SideIterator sideIterEnd = m_uInfo->template end<side_type>(si);
613 for( ;sideIter !=sideIterEnd; sideIter++)
614 {
615 side_type* side = *sideIter;
616 if (pbm && pbm->is_slave(side)) continue;
617 dimMat Tij;
618 MathVector<dim> uValue;
619 if (useGridFunction)
620 for (int d=0;d<dim;d++){
621 u->dof_indices(side, d, multInd);
622 uValue[d]=DoFRef(*u,multInd[0]);
623 }
624 else
625 for (int d=0;d<dim;d++) uValue[d] = (*aaU)[side][d];
626 size_t ind=0;
627 for (int d1=0;d1 < dim;d1++)
628 for (int d2=d1;d2 < dim;d2++){
629 Tij[ind] = uValue[d1]*uValue[d2];
630 ind++;
631 }
632 Tij*=factor;
633 aaResult[side]+=Tij;
634 }
635 }
636}
637
638template<typename TGridFunction>
640 // get domain of grid function
641 domain_type& domain = *m_u->domain().get();
642 SetAttachmentValues(m_acTurbulentViscosity, m_grid->template begin<side_type>(), m_grid->template end<side_type>(), 0);
643
644 // coord and vertex array
645 // MathVector<dim> coCoord[domain_traits<dim>::MaxNumVerticesOfElem];
646 // Vertex* vVrt[domain_traits<dim>::MaxNumVerticesOfElem];
647
648 // assemble deformation tensor fluxes
649 this->assembleDeformationTensor(m_acDeformation,m_acVolume,m_u,NULL);
650 // compute turbulent viscosity , loop over sides
651 for(int si = 0; si < domain.subset_handler()->num_subsets(); ++si)
652 {
653 if ((m_turbZeroSg.size()!=0) && (m_turbZeroSg.contains(si)==true)) continue;
654 SideIterator sideIter = m_u->template begin<side_type>(si);
655 SideIterator sideIterEnd = m_u->template end<side_type>(si);
656 for( ;sideIter !=sideIterEnd; sideIter++)
657 {
658 // get Elem
659 side_type* side = *sideIter;
660 if (m_pbm && m_pbm->is_slave(side)) continue;
661 number delta = m_acVolume[side];
662 // for possible other choices of delta see Fr�hlich p 160
663 delta = pow(delta,(number)1.0/(number)dim);
664 m_acTurbulentViscosity[side] = m_c*m_c * delta*delta * m_acDeformation[side].fnorm();
665 }
666 }
667 if (m_bAdaptive){
668 if (dim==2) constrainingSideAveraging<TGridFunction,side_type,ConstrainingEdge,number>(m_acTurbulentViscosity,m_u);
669 else {
670 constrainingSideAveraging<TGridFunction,side_type,ConstrainingTriangle,number>(m_acTurbulentViscosity,m_u);
671 constrainingSideAveraging<TGridFunction,side_type,ConstrainingQuadrilateral,number>(m_acTurbulentViscosity,m_u);
672 }
673 }
674 // transfer attachment data to lower levels
675 this->transferToLowerLevels(m_acTurbulentViscosity,*m_spApproxSpace);
676}
677
678template<typename TGridFunction>
680 // get domain of grid function
681 domain_type& domain = *m_u->domain().get();
682
683 // get position accessor
684 // for debug
685 typedef typename domain_type::position_accessor_type position_accessor_type;
686 // for debug
687 const position_accessor_type& posAcc = domain.position_accessor();
688
689 // initialize attachment values with 0
690// SetAttachmentValues(m_acDeformation , m_grid->template begin<side_type>(), m_grid->template end<side_type>(), 0);
691 SetAttachmentValues(m_acTurbulentViscosity, m_grid->template begin<side_type>(), m_grid->template end<side_type>(), 0);
692// SetAttachmentValues(m_acVolume,m_grid->template begin<side_type>(), m_grid->template end<side_type>(), 0);
693 SetAttachmentValues(m_acTurbulentCNew,m_grid->template begin<side_type>(), m_grid->template end<side_type>(), 0);
694// SetAttachmentValues(m_acVolumeHat,m_grid->template begin<side_type>(), m_grid->template end<side_type>(), 0);
695// SetAttachmentValues(m_acUHat,m_grid->template begin<side_type>(), m_grid->template end<side_type>(), 0);
696// SetAttachmentValues(m_acDeformationHat,m_grid->template begin<side_type>(), m_grid->template end<side_type>(), 0);
697// SetAttachmentValues(m_acLij,m_grid->template begin<side_type>(), m_grid->template end<side_type>(), 0);
698 SetAttachmentValues(m_acMij,m_grid->template begin<side_type>(), m_grid->template end<side_type>(), 0);
699
700 // compute Lij term \hat{u_i u_j} - \hat{u_i} \hat{u_j}
701 // \hat{u}
702 this->elementFilter(m_acUHat,m_acVolumeHat,m_u,&m_acUHat);
703 // use Mij attachment to store first Lij part
704 // u_i u_j
705 this->addUiUjTerm(m_acMij,1.0,m_u,NULL);
706 // \hat{u_i u_j}
707 this->elementFilter(m_acLij,m_acVolumeHat,SPNULL,&m_acMij);
708 // \hat{u_i u_j} - \hat{u_i} \hat{u_j}
709 this->addUiUjTerm(m_acLij,-1.0,SPNULL,&m_acUHat);
710
711 // Mij term
712 // first term |\hat{S}| \hat{S}
713 // assemble \hat{S} using \hat{u}
714 this->assembleDeformationTensor(m_acDeformationHat,m_acVolume,SPNULL,&m_acUHat);
715 // normalize \hat{S}
716 this->scaleTensorByNorm(m_acDeformationHat);
717 // Mij second term \hat{|S|S}
718 // compute S
719 this->assembleDeformationTensor(m_acDeformation,m_acVolume,m_u,NULL);
720 // compute |S| S
721 this->scaleTensorByNorm(m_acDeformation);
722 // filter |S| S
723 //for debug UG_LOG("------------------------------------------------------\n");
724 this->elementFilter(m_acMij,m_acVolumeHat,SPNULL,&m_acDeformation);
725
726 // complete Mij term computation by scaling and adding the two terms,
727 // solve the local least squares problem and compute local c and local turbulent viscosity
728 for(int si = 0; si < domain.subset_handler()->num_subsets(); ++si)
729 {
730 if ((m_turbZeroSg.size()!=0) && (m_turbZeroSg.contains(si)==true)) continue;
731 SideIterator sideIter = m_u->template begin<side_type>(si);
732 SideIterator sideIterEnd = m_u->template end<side_type>(si);
733 for( ;sideIter !=sideIterEnd; sideIter++){
734 side_type* side = *sideIter;
735 if (m_pbm && m_pbm->is_slave(side)) continue;
736 // use fixed ratio \hat{delta} / delta or compute it from filter width
737 number kappa;
738 if (m_fixedRatio) kappa = m_kappa;
739 else {
740 number delta = m_acVolume[side];
741 delta = pow(delta,(number)1.0/(number)dim);
742 number deltaHat = m_acVolumeHat[side];
743 deltaHat = pow(deltaHat,(number)1.0/(number)dim);
744 kappa = deltaHat /(number) delta;
745 }
746 m_acDeformationHat[side] *= -kappa*kappa;
747 m_acMij[side] += m_acDeformationHat[side];
748 // compute local c
749 // solve least squares problem
750 number c = 0;
751 for (int d1=0;d1<dim;d1++)
752 for (int d2=0;d2<dim;d2++)
753 c += m_acLij[side](d1,d2)*m_acMij[side](d1,d2);
754 number denom=0;
755 //for debug UG_LOG("c=" << c << "\n");
756 for (int d1=0;d1<dim;d1++)
757 for (int d2=0;d2<dim;d2++)
758 denom += m_acMij[side](d1,d2)*m_acMij[side](d1,d2);
759 if (denom>1e-15)
760 c/=(number)(2.0*denom);
761 else c=0;
762
763// UG_LOG("c=" << c << "\n");
764
765 if (m_spaceFilter==false){
766 if (m_timeFilter==false){
767 m_acTurbulentViscosity[side] = c * m_acDeformation[side].fnorm();
768 UG_LOG("nu_t = " << m_acTurbulentViscosity[side] << " c = " << c << " kappa = " << kappa << " co=[" << 0.5*(posAcc[side->vertex(0)][0] + posAcc[side->vertex(1)][0]) << "," << 0.5*(posAcc[side->vertex(0)][1] + posAcc[side->vertex(1)][1]) << "]\n");
769 } else {
770 m_acTurbulentC[side]= (m_timeFilterEps * c + (1-m_timeFilterEps)*m_acTurbulentC[side]);
771 m_acTurbulentViscosity[side] = m_acTurbulentC[side] * m_acDeformation[side].fnorm();
772 }
773 if (m_acTurbulentViscosity[side]+m_viscosityNumber<m_small) m_acTurbulentViscosity[side] = m_viscosityNumber + m_small; }
774 else{
775 // store c in viscosity attachment
776 m_acTurbulentViscosity[side] = c;
777 }
778 }
779 }
780 if (m_spaceFilter==true){
781 // filter c
782 if (m_timeFilter==false)
783 // c has been stored in viscosity attachment
784 this->elementFilter(m_acTurbulentC,m_acVolumeHat,SPNULL,&m_acTurbulentViscosity);
785 else
786 this->elementFilter(m_acTurbulentCNew,m_acVolumeHat,SPNULL,&m_acTurbulentViscosity);
787 // compute turbulent viscosity
788 for(int si = 0; si < domain.subset_handler()->num_subsets(); ++si)
789 {
790 //for debug UG_LOG("si = " << si << "\n");
791 if ((m_turbZeroSg.size()!=0) && (m_turbZeroSg.contains(si)==true)) continue;
792 SideIterator sideIter = m_u->template begin<side_type>(si);
793 SideIterator sideIterEnd = m_u->template end<side_type>(si);
794 for( ;sideIter !=sideIterEnd; sideIter++){
795 side_type* side = *sideIter;
796 if (m_pbm && m_pbm->is_slave(side)) continue;
797 if (m_timeFilter==true)
798 // time averaging, note that c has been stored in m_acVolumeHat
799 m_acTurbulentC[side]= (m_timeFilterEps * m_acTurbulentCNew[side] + (1-m_timeFilterEps)*m_acTurbulentC[side]);
800 m_acTurbulentViscosity[side] = m_acTurbulentC[side] * m_acDeformation[side].fnorm();
801 UG_LOG("nu_t = " << m_acTurbulentViscosity[side] << " c = " << m_acTurbulentC[side] << " co=[" << 0.5*(posAcc[side->vertex(0)][0] + posAcc[side->vertex(1)][0]) << "," << 0.5*(posAcc[side->vertex(0)][1] + posAcc[side->vertex(1)][1]) << "]\n");
802 if (m_acTurbulentViscosity[side]+m_viscosityNumber<m_small) m_acTurbulentViscosity[side] = m_viscosityNumber+m_small;
803 // for debug
804 // UG_LOG("nu_t = " << m_acTurbulentViscosity[side] << " c = " << m_acTurbulentC[side] << " co=[" << 0.5*(posAcc[side->vertex(0)][0] + posAcc[side->vertex(1)][0]) << "," << 0.5*(posAcc[side->vertex(0)][1] + posAcc[side->vertex(1)][1]) << "]\n");
805 }
806 }
807 }
808 if (m_bAdaptive){
809 if (dim==2) constrainingSideAveraging<TGridFunction,side_type,ConstrainingEdge,number>(m_acTurbulentViscosity,m_u);
810 else {
811 constrainingSideAveraging<TGridFunction,side_type,ConstrainingTriangle,number>(m_acTurbulentViscosity,m_u);
812 constrainingSideAveraging<TGridFunction,side_type,ConstrainingQuadrilateral,number>(m_acTurbulentViscosity,m_u);
813 }
814 }
815 // transfer attachment data to lower levels
816 this->transferToLowerLevels(m_acTurbulentViscosity,*m_spApproxSpace);
817}
818
819} // namespace NavierStokes
820} // end namespace ug
821
822#endif /* __H__UG__NAVIER_STOKES_TURBULENT_VISCOSITY_DATA_IMPL_H_ */
parameterString s
Definition Biogas.lua:2
function NavierStokes(fcts, subsets, discType)
T * get()
SmartPtr< TDomain > domain()
size_t node_id() const
const MathVector< worldDim > & normal() const
const MathVector< worldDim > & normal() const
number shape(size_t sh) const
size_t num_corners() const
size_t node_id() const
const MathVector< worldDim > & global_corner(size_t co) const
const std::vector< BF > & bf(int subsetIndex) const
size_t num_scvf() const
size_t num_scv() const
const SCV & scv(size_t i) const
void add_boundary_subset(int subsetIndex)
size_t num_sh() const
void update(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
size_t num_bf() const
const SCVF & scvf(size_t i) const
size_t num_constrained_dofs() const
void update_hanging(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL, bool keepSCV=false, bool keepSCVF=false)
virtual void global_to_local(MathVector< dim > &locPos, const MathVector< worldDim > &globPos, const size_t maxIter=1000, const number tol=1e-10) const=0
SmartPtr< MultiGrid > multi_grid()
size_t num_levels() const
SmartPtr< DoFDistribution > dof_distribution(const GridLevel &gl, bool bCreate=true)
virtual void shapes(std::vector< std::vector< shape_type > > &vvShape, const std::vector< MathVector< dim > > &vLocPos) const=0
TGridFunction::domain_type domain_type
domain type
Definition turbulent_viscosity_fvcr.h:418
domain_type::position_accessor_type position_accessor_type
position accessor type
Definition turbulent_viscosity_fvcr.h:424
elem_type::side side_type
side type
Definition turbulent_viscosity_fvcr.h:436
TGridFunction::template traits< side_type >::const_iterator SideIterator
side iterator
Definition turbulent_viscosity_fvcr.h:442
void update()
Definition turbulent_viscosity_fvcr_impl.h:679
TGridFunction::domain_type domain_type
domain type
Definition turbulent_viscosity_fvcr.h:301
TGridFunction::template traits< side_type >::const_iterator SideIterator
side iterator
Definition turbulent_viscosity_fvcr.h:325
void update()
Definition turbulent_viscosity_fvcr_impl.h:639
elem_type::side side_type
side type
Definition turbulent_viscosity_fvcr.h:319
TGridFunction::template traits< side_type >::const_iterator SideIterator
side iterator
Definition turbulent_viscosity_fvcr.h:86
TGridFunction::domain_type domain_type
domain type
Definition turbulent_viscosity_fvcr.h:65
void transferToLowerLevels(aSideNumber &aaData, ApproximationSpace< domain_type > &approximationSpace)
Definition turbulent_viscosity_fvcr_impl.h:40
TGridFunction::template dim_traits< dim >::const_iterator ElemIterator
element iterator
Definition turbulent_viscosity_fvcr.h:83
void addUiUjTerm(aSideTensor &aaDefTensor, const number factor, SmartPtr< TGridFunction > u, aSideDimVector *aaU)
Definition turbulent_viscosity_fvcr_impl.h:598
elem_type::side side_type
side type
Definition turbulent_viscosity_fvcr.h:80
void elementFilter(PeriodicAttachmentAccessor< side_type, Attachment< VType > > &aaUHat, aSideNumber &aaVol, SmartPtr< TGridFunction > u, PeriodicAttachmentAccessor< side_type, Attachment< VType > > *aaU)
Definition turbulent_viscosity_fvcr_impl.h:109
void scaleTensorByNorm(aSideTensor &aaTensor)
Definition turbulent_viscosity_fvcr_impl.h:575
void fillAttachment(aSideDimVector &aaU, SmartPtr< TGridFunction > u)
Definition turbulent_viscosity_fvcr_impl.h:66
void scvFilter(PeriodicAttachmentAccessor< side_type, Attachment< VType > > &aaUHat, aSideNumber &aaVol, SmartPtr< TGridFunction > u, PeriodicAttachmentAccessor< side_type, Attachment< VType > > *aaU)
Definition turbulent_viscosity_fvcr_impl.h:262
void assembleDeformationTensor(aSideTensor &aaDefTensor, aSideNumber &aaVol, SmartPtr< TGridFunction > u, aSideDimVector *aaU)
Definition turbulent_viscosity_fvcr_impl.h:406
domain_type::position_accessor_type position_accessor_type
position accessor type
Definition turbulent_viscosity_fvcr.h:71
TGridFunction::template dim_traits< dim >::grid_base_object elem_type
element type
Definition turbulent_viscosity_fvcr.h:77
bool is_slave(TElem *) const
SmartPtr< TSubsetHandler > subset_handler()
SmartPtr< TGrid > grid()
size_t num_subsets() const
void SetAttachmentValues(TAttachmentAccessor &aaVal, TIter elemsBegin, TIter elemsEnd, const TVal &val)
const NullSmartPtr SPNULL
#define UG_ASSERT(expr, msg)
#define UG_LOG(msg)
double number
void assignValue(number &v, size_t ind, number value)
Definition turbulent_viscosity_fvcr_impl.h:85
ReferenceObjectID
const number & DoFRef(const TMatrix &mat, const DoFIndex &iInd, const DoFIndex &jInd)