My Project
programmer's documentation
Loading...
Searching...
No Matches
cs_hho_builder.h
Go to the documentation of this file.
1#ifndef __CS_HHO_BUILDER_H__
2#define __CS_HHO_BUILDER_H__
3
4/*============================================================================
5 * Low-level functions and structures used to build the algebraic system with
6 * a cellwise process when Hybrid High Order schemes are set for the space
7 * discretization
8 *============================================================================*/
9
10/*
11 This file is part of Code_Saturne, a general-purpose CFD tool.
12
13 Copyright (C) 1998-2019 EDF S.A.
14
15 This program is free software; you can redistribute it and/or modify it under
16 the terms of the GNU General Public License as published by the Free Software
17 Foundation; either version 2 of the License, or (at your option) any later
18 version.
19
20 This program is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
22 FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
23 details.
24
25 You should have received a copy of the GNU General Public License along with
26 this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
27 Street, Fifth Floor, Boston, MA 02110-1301, USA.
28*/
29
30/*----------------------------------------------------------------------------
31 * Local headers
32 *----------------------------------------------------------------------------*/
33
34#include "cs_base.h"
35#include "cs_basis_func.h"
36#include "cs_cdo_connect.h"
37#include "cs_sdm.h"
38#include "cs_xdef.h"
39
40/*----------------------------------------------------------------------------*/
41
43
44/*============================================================================
45 * Macro definitions
46 *============================================================================*/
47
48/*============================================================================
49 * Type definitions
50 *============================================================================*/
51
52/* Cellwise builder for HHO discretization */
53typedef struct {
54
55 /* Current and maximal number of face basis allocated. This number is equal
56 to the max. number of faces for a cell */
57 short int n_face_basis;
59
60 cs_basis_func_t **face_basis; /* P_(d-1)^k polynomial basis */
61 cs_basis_func_t *cell_basis; /* P_d^k polynomial basis */
62 cs_basis_func_t *grad_basis; /* P_d^(k+1) \ P_d^0 polynomial basis */
63
64 cs_sdm_t *grad_reco_op; /* Gradient operator; Rectangular matrix */
65
66 /* Temporary matrices defined by blocks */
67 cs_sdm_t *tmp; /* Temporary block matrix (fs x ts) */
68 cs_sdm_t *bf_t; /* Transposed of Bf (used in stabilization) */
69 cs_sdm_t *jstab; /* Stabilization part related to a face */
70
72
73/*============================================================================
74 * Public function prototypes
75 *============================================================================*/
76
77/*----------------------------------------------------------------------------*/
86/*----------------------------------------------------------------------------*/
87
89cs_hho_builder_create(int order,
90 int n_fc);
91
92/*----------------------------------------------------------------------------*/
98/*----------------------------------------------------------------------------*/
99
100void
102
103/*----------------------------------------------------------------------------*/
114/*----------------------------------------------------------------------------*/
115
116void
119 cs_hho_builder_t *hhob);
120
121/*----------------------------------------------------------------------------*/
129/*----------------------------------------------------------------------------*/
130
131static inline void
134 cs_hho_builder_t *hhob)
135{
136 if (hhob == NULL)
137 return;
138
139 /* Sanity checks */
140 assert(cm != NULL);
141
142 /* Setup cell basis functions */
143 hhob->cell_basis->setup(hhob->cell_basis, cm, 0, cm->xc, cb);
144 hhob->n_face_basis = 0;
145}
146
147/*----------------------------------------------------------------------------*/
159/*----------------------------------------------------------------------------*/
160
161void
164 cs_hho_builder_t *hhob);
165
166/*----------------------------------------------------------------------------*/
175/*----------------------------------------------------------------------------*/
176
177void
180 cs_hho_builder_t *hhob);
181
182/*----------------------------------------------------------------------------*/
196/*----------------------------------------------------------------------------*/
197
198void
200 const cs_cell_mesh_t *cm,
201 cs_real_t t_eval,
203 cs_hho_builder_t *hhob,
204 cs_real_t red[]);
205
206/*----------------------------------------------------------------------------*/
222/*----------------------------------------------------------------------------*/
223
224void
226 const cs_cell_mesh_t *cm,
227 cs_real_t t_eval,
229 cs_hho_builder_t *hhob,
230 cs_real_t red[]);
231
232/*----------------------------------------------------------------------------*/
245/*----------------------------------------------------------------------------*/
246
247void
249 short int f,
250 const cs_cell_mesh_t *cm,
251 cs_real_t t_eval,
253 cs_hho_builder_t *hhob,
254 cs_real_t res[]);
255
256/*----------------------------------------------------------------------------*/
269/*----------------------------------------------------------------------------*/
270
271void
273 short int f,
274 const cs_cell_mesh_t *cm,
275 cs_real_t t_eval,
277 cs_hho_builder_t *hhob,
278 cs_real_t res[]);
279
280/*----------------------------------------------------------------------------*/
281
283
284#endif /* __CS_HHO_BUILDER_H__ */
#define BEGIN_C_DECLS
Definition cs_defs.h:467
double cs_real_t
Floating-point value.
Definition cs_defs.h:302
#define END_C_DECLS
Definition cs_defs.h:468
void cs_hho_builder_compute_grad_reco(const cs_cell_mesh_t *cm, cs_cell_builder_t *cb, cs_hho_builder_t *hhob)
Compute the gradient operator stemming from the relation stiffness * grad_op = rhs where stiffness is...
Definition cs_hho_builder.c:892
void cs_hho_builder_free(cs_hho_builder_t **p_builder)
Free a cs_hho_builder_t structure.
Definition cs_hho_builder.c:802
void cs_hho_builder_compute_dirichlet_v(const cs_xdef_t *def, short int f, const cs_cell_mesh_t *cm, cs_real_t t_eval, cs_cell_builder_t *cb, cs_hho_builder_t *hhob, cs_real_t res[])
Compute the projection of the Dirichlet boundary conditions onto the polynomial spaces on faces....
Definition cs_hho_builder.c:1940
static void cs_hho_builder_cellbasis_setup(const cs_cell_mesh_t *cm, cs_cell_builder_t *cb, cs_hho_builder_t *hhob)
Set-up the basis functions related to a cell only.
Definition cs_hho_builder.h:132
void cs_hho_builder_compute_dirichlet(const cs_xdef_t *def, short int f, const cs_cell_mesh_t *cm, cs_real_t t_eval, cs_cell_builder_t *cb, cs_hho_builder_t *hhob, cs_real_t res[])
Compute the projection of the Dirichlet boundary conditions onto the polynomial spaces on faces.
Definition cs_hho_builder.c:1813
cs_hho_builder_t * cs_hho_builder_create(int order, int n_fc)
Allocate a cs_hho_builder_t structure.
Definition cs_hho_builder.c:749
void cs_hho_builder_reduction_from_analytic(const cs_xdef_t *def, const cs_cell_mesh_t *cm, cs_real_t t_eval, cs_cell_builder_t *cb, cs_hho_builder_t *hhob, cs_real_t red[])
Compute the reduction onto the polynomial spaces (cell and faces) of a function defined by an analyti...
Definition cs_hho_builder.c:1435
void cs_hho_builder_diffusion(const cs_cell_mesh_t *cm, cs_cell_builder_t *cb, cs_hho_builder_t *hhob)
Compute the diffusion operator. The gradient reconstruction operator has to be built just before this...
Definition cs_hho_builder.c:1128
void cs_hho_builder_reduction_from_analytic_v(const cs_xdef_t *def, const cs_cell_mesh_t *cm, cs_real_t t_eval, cs_cell_builder_t *cb, cs_hho_builder_t *hhob, cs_real_t red[])
Compute the reduction onto the polynomial spaces (cell and faces) of a function defined by an analyti...
Definition cs_hho_builder.c:1609
void cs_hho_builder_cellwise_setup(const cs_cell_mesh_t *cm, cs_cell_builder_t *cb, cs_hho_builder_t *hhob)
Set-up the basis functions related to a cell, its gradient and to the faces of this cell....
Definition cs_hho_builder.c:842
Definition cs_basis_func.h:181
cs_basis_func_setup_t * setup
Definition cs_basis_func.h:199
Set of local and temporary buffers useful for building the algebraic system with a cellwise process....
Definition cs_cdo_local.h:56
Set of local quantities and connectivities related to a mesh cell This is a key structure for all cel...
Definition cs_cdo_local.h:146
cs_real_3_t xc
Definition cs_cdo_local.h:158
Definition cs_hho_builder.h:53
cs_sdm_t * tmp
Definition cs_hho_builder.h:67
cs_basis_func_t ** face_basis
Definition cs_hho_builder.h:60
cs_basis_func_t * grad_basis
Definition cs_hho_builder.h:62
short int n_max_face_basis
Definition cs_hho_builder.h:58
short int n_face_basis
Definition cs_hho_builder.h:57
cs_basis_func_t * cell_basis
Definition cs_hho_builder.h:61
cs_sdm_t * grad_reco_op
Definition cs_hho_builder.h:64
cs_sdm_t * jstab
Definition cs_hho_builder.h:69
cs_sdm_t * bf_t
Definition cs_hho_builder.h:68
Structure storing medata for defining a quantity in a very flexible way.
Definition cs_xdef.h:126