OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
mesh_worker_workers.h
Go to the documentation of this file.
1 // ------------------------------------------------------------------
2 //
3 // $Id: mesh_worker_workers.h 1348 2013-08-16 00:45:07Z secanell $
4 //
5 // Copyright (C) 2006, 2007, 2008, 2009 by Guido Kanschat.
6 // Copyright (C) 2012 by Valentin N. Zingan, University of Alberta.
7 //
8 // This file is subject to QPL and may not be distributed
9 // without copyright and license information. Please refer
10 // to the file deal.II/doc/license.html for the text and
11 // further information on this license.
12 //
13 // ------------------------------------------------------------------
14 
15 #ifndef _DEALII_APPFRAME_MESH_WORKER_WORKERS_H_
16 #define _DEALII_APPFRAME_MESH_WORKER_WORKERS_H_
17 
19 
20 // ################################################################################################
21 
36 // ##### CATEGORY-2 ################################################# WORKER_OBJECTS ##############
37 //
38 // ##### MeshWorker::WorkerObjects::LocalWorker #####
39 // ##### MeshWorker::WorkerObjects::IntegrationWorker #####
40 
41 namespace MeshWorker // # 2
42 {
43 
44 namespace WorkerObjects
45 {
46 
47 // --- 1 ---
48 
65 template<int dim>
67 {
68 public:
69 
79  {
80  boundary_fluxes = true;
81  interior_fluxes = true;
82  }
83 
87  void cell(const InfoObjects::DoFInfo<dim>& cell);
88 
92  void bdry(const InfoObjects::DoFInfo<dim>& face);
93 
97  void face(const InfoObjects::DoFInfo<dim>& face1,
98  const InfoObjects::DoFInfo<dim>& face2);
99 
105 
111 };
112 
113 // --- 2 ---
114 
121 template<int dim>
122 class IntegrationWorker : public LocalWorker<dim>
123 {
124 public:
125 
132 
151  void initialize_selectors(const VectorSelector& cs,
152  const VectorSelector& bs,
153  const VectorSelector& fs);
154 
160  void add_selector(const std::string& name,
161  bool values,
162  bool gradients,
163  bool hessians,
164  bool cell,
165  bool bdry,
166  bool face);
167 
171  void add_update_flags(const UpdateFlags flags,
172  bool cell = true,
173  bool bdry = true,
174  bool face = true,
175  bool ngbr = true)
176  {
177  if(cell) cell_flags |= flags;
178  if(bdry) bdry_flags |= flags;
179  if(face) face_flags |= flags;
180  if(ngbr) ngbr_flags |= flags;
181  }
182 
192  void initialize_gauss_quadrature(unsigned int n_cell_points,
193  unsigned int n_bdry_points,
194  unsigned int n_face_points)
195  {
196  cell_quadrature = QGauss<dim>(n_cell_points);
197  bdry_quadrature = QGauss<dim-1>(n_bdry_points);
198  face_quadrature = QGauss<dim-1>(n_face_points);
199  }
200 
208 
216 
224 
229  Quadrature<dim> cell_quadrature;
230 
235  Quadrature<dim-1> bdry_quadrature;
236 
241  Quadrature<dim-1> face_quadrature;
242 
250  UpdateFlags cell_flags;
251 
261  UpdateFlags bdry_flags;
262 
272  UpdateFlags face_flags;
273 
284  UpdateFlags ngbr_flags;
285 };
286 
287 // --- implementation ---
288 
289 template<int dim>
290 inline
292 :
293 LocalWorker<dim>()
294 {
295  cell_flags = update_JxW_values;
296  bdry_flags = UpdateFlags(update_normal_vectors | update_JxW_values);
298  ngbr_flags = update_default;
299 }
300 
301 template<int dim>
302 inline
303 void
305  const VectorSelector& bs,
306  const VectorSelector& fs)
307 {
308  cell_selector = cs;
309  bdry_selector = bs;
310  face_selector = fs;
311 
312  if(cell_selector.has_values() != 0) cell_flags |= update_values;
313  if(cell_selector.has_gradients() != 0) cell_flags |= update_gradients;
314  if(cell_selector.has_hessians() != 0) cell_flags |= update_hessians;
315 
316  if(bdry_selector.has_values() != 0) bdry_flags |= update_values;
317  if(bdry_selector.has_gradients() != 0) bdry_flags |= update_gradients;
318  if(bdry_selector.has_hessians() != 0) bdry_flags |= update_hessians;
319 
320  if(face_selector.has_values() != 0) face_flags |= update_values;
321  if(face_selector.has_gradients() != 0) face_flags |= update_gradients;
322  if(face_selector.has_hessians() != 0) face_flags |= update_hessians;
323 
324  if(face_selector.has_values() != 0) ngbr_flags |= update_values;
325  if(face_selector.has_gradients() != 0) ngbr_flags |= update_gradients;
326  if(face_selector.has_hessians() != 0) ngbr_flags |= update_hessians;
327 }
328 
329 template<int dim>
330 inline
331 void
332 IntegrationWorker<dim>::add_selector(const std::string& name,
333  bool values,
334  bool gradients,
335  bool hessians,
336  bool cell,
337  bool bdry,
338  bool face)
339 {
340  if(cell) cell_selector.add_vector(name, values, gradients, hessians);
341  if(bdry) bdry_selector.add_vector(name, values, gradients, hessians);
342  if(face) face_selector.add_vector(name, values, gradients, hessians);
343 
344  if(cell_selector.has_values() != 0) cell_flags |= update_values;
345  if(cell_selector.has_gradients() != 0) cell_flags |= update_gradients;
346  if(cell_selector.has_hessians() != 0) cell_flags |= update_hessians;
347 
348  if(bdry_selector.has_values() != 0) bdry_flags |= update_values;
349  if(bdry_selector.has_gradients() != 0) bdry_flags |= update_gradients;
350  if(bdry_selector.has_hessians() != 0) bdry_flags |= update_hessians;
351 
352  if(face_selector.has_values() != 0) face_flags |= update_values;
353  if(face_selector.has_gradients() != 0) face_flags |= update_gradients;
354  if(face_selector.has_hessians() != 0) face_flags |= update_hessians;
355 
356  if(face_selector.has_values() != 0) ngbr_flags |= update_values;
357  if(face_selector.has_gradients() != 0) ngbr_flags |= update_gradients;
358  if(face_selector.has_hessians() != 0) ngbr_flags |= update_hessians;
359 }
360 
361 } // namespace WorkerObjects
362 
363 } // namespace MeshWorker # 2
364 
365 #endif