My Project
programmer's documentation
Loading...
Searching...
No Matches
cs_hodge.h
Go to the documentation of this file.
1#ifndef __CS_HODGE_H__
2#define __CS_HODGE_H__
3
4/*============================================================================
5 * Build discrete Hodge operators
6 *============================================================================*/
7
8/*
9 This file is part of Code_Saturne, a general-purpose CFD tool.
10
11 Copyright (C) 1998-2019 EDF S.A.
12
13 This program is free software; you can redistribute it and/or modify it under
14 the terms of the GNU General Public License as published by the Free Software
15 Foundation; either version 2 of the License, or (at your option) any later
16 version.
17
18 This program is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
20 FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
21 details.
22
23 You should have received a copy of the GNU General Public License along with
24 this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
25 Street, Fifth Floor, Boston, MA 02110-1301, USA.
26*/
27
28/*----------------------------------------------------------------------------
29 * Standard C library headers
30 *----------------------------------------------------------------------------*/
31
32/*----------------------------------------------------------------------------
33 * Local headers
34 *----------------------------------------------------------------------------*/
35
36#include "cs_cdo_connect.h"
37#include "cs_cdo_local.h"
38#include "cs_cdo_quantities.h"
39#include "cs_param.h"
40#include "cs_param_cdo.h"
41#include "cs_property.h"
42
43/*----------------------------------------------------------------------------*/
44
46
47/*============================================================================
48 * Macro and type definitions
49 *============================================================================*/
50
51/*----------------------------------------------------------------------------*/
60/*----------------------------------------------------------------------------*/
61
62typedef void
63(cs_hodge_t)(const cs_param_hodge_t h_info,
64 const cs_cell_mesh_t *cm,
66
67/*============================================================================
68 * Public function definitions
69 *============================================================================*/
70
71/*----------------------------------------------------------------------------*/
81/*----------------------------------------------------------------------------*/
82
83void
85 const cs_cell_mesh_t *cm,
87
88/*----------------------------------------------------------------------------*/
98/*----------------------------------------------------------------------------*/
99
100void
102 const cs_cell_mesh_t *cm,
104
105/*----------------------------------------------------------------------------*/
115/*----------------------------------------------------------------------------*/
116
117void
119 const cs_cell_mesh_t *cm,
121
122/*----------------------------------------------------------------------------*/
131/*----------------------------------------------------------------------------*/
132
133void
135 const cs_cell_mesh_t *cm,
137
138/*----------------------------------------------------------------------------*/
147/*----------------------------------------------------------------------------*/
148
149void
151 const cs_cell_mesh_t *cm,
153
154/*----------------------------------------------------------------------------*/
164/*----------------------------------------------------------------------------*/
165
166void
168 const cs_cell_mesh_t *cm,
170
171/*----------------------------------------------------------------------------*/
180/*----------------------------------------------------------------------------*/
181
182void
184 const cs_cell_mesh_t *cm,
186
187/*----------------------------------------------------------------------------*/
197/*----------------------------------------------------------------------------*/
198
199void
201 const cs_cell_mesh_t *cm,
203
204/*----------------------------------------------------------------------------*/
214/*----------------------------------------------------------------------------*/
215
216void
218 const cs_cell_mesh_t *cm,
220
221/*----------------------------------------------------------------------------*/
230/*----------------------------------------------------------------------------*/
231
232void
234 const cs_cell_mesh_t *cm,
236
237/*----------------------------------------------------------------------------*/
247/*----------------------------------------------------------------------------*/
248
249void
251 const cs_cell_mesh_t *cm,
253
254/*----------------------------------------------------------------------------*/
264/*----------------------------------------------------------------------------*/
265
266void
268 const cs_cell_mesh_t *cm,
270
271/*----------------------------------------------------------------------------*/
281/*----------------------------------------------------------------------------*/
282
283void
285 const cs_cell_mesh_t *cm,
287
288/*----------------------------------------------------------------------------*/
298/*----------------------------------------------------------------------------*/
299
300void
302 const cs_cell_mesh_t *cm,
304
305/*----------------------------------------------------------------------------*/
315/*----------------------------------------------------------------------------*/
316
317void
319 const cs_cell_mesh_t *cm,
321
322/*----------------------------------------------------------------------------*/
332/*----------------------------------------------------------------------------*/
333
334void
336 const cs_cell_mesh_t *cm,
338
339/*----------------------------------------------------------------------------*/
349/*----------------------------------------------------------------------------*/
350
351void
353 const cs_cell_mesh_t *cm,
355
356/*----------------------------------------------------------------------------*/
366/*----------------------------------------------------------------------------*/
367
368void
370 const cs_cell_mesh_t *cm,
372
373/*----------------------------------------------------------------------------*/
383/*----------------------------------------------------------------------------*/
384
385void
387 const cs_cell_mesh_t *cm,
389
390/*----------------------------------------------------------------------------*/
403/*----------------------------------------------------------------------------*/
404
405void
406cs_hodge_matvec(const cs_cdo_connect_t *connect,
407 const cs_cdo_quantities_t *quant,
408 const cs_param_hodge_t h_info,
409 const cs_property_t *pty,
410 const double in_vals[],
411 cs_real_t t_eval,
412 double result[]);
413
414/*----------------------------------------------------------------------------*/
423/*----------------------------------------------------------------------------*/
424
425void
427 cs_sdm_t *hf);
428
429/*----------------------------------------------------------------------------*/
430
432
433#endif /* __CS_HODGE_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
@ fm
Definition cs_field_pointer.h:134
void cs_hodge_edfp_cost_get_opt(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local Hodge operator for a given cell using the COST algo. Hodge op. from dual edges to prima...
Definition cs_hodge.c:2517
void cs_hodge_compute_wbs_surfacic(const cs_face_mesh_t *fm, cs_sdm_t *hf)
Compute the hodge operator related to a face (i.e. a mass matrix with unity property) using a Whitney...
Definition cs_hodge.c:2854
void cs_hodge_fb_voro_get_stiffness(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local stiffness matrix using the Voronoi algorithm The computed matrix is stored in cb->loc C...
Definition cs_hodge.c:753
void cs_hodge_fped_voro_get(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local Hodge operator for a given cell using VORONOI algo. Hodge op. from primal faces to dual...
Definition cs_hodge.c:2247
void cs_hodge_epfd_voro_get(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local Hodge operator for a given cell using VORONOI algo. Hodge op. from primal edges to dual...
Definition cs_hodge.c:2105
void cs_hodge_fb_get_mass(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local Hodge operator on a given cell which is equivalent of a mass matrix....
Definition cs_hodge.c:1756
void cs_hodge_edfp_voro_get(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local Hodge operator for a given cell using VORONOI algo. Hodge op. from dual edges to primal...
Definition cs_hodge.c:2382
void cs_hodge_vb_cost_get_stiffness(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local stiffness matrix using the generic COST algo. The computed matrix is stored in cb->loc ...
Definition cs_hodge.c:1164
void cs_hodge_vb_cost_get_iso_stiffness(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local stiffness matrix using the generic COST algo. Case of CDO vertex-based schemes and isot...
Definition cs_hodge.c:814
void cs_hodge_matvec(const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_param_hodge_t h_info, const cs_property_t *pty, const double in_vals[], cs_real_t t_eval, double result[])
Compute cellwise a discrete hodge operator and multiple it with a vector.
Definition cs_hodge.c:2598
void cs_hodge_fped_cost_get(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local Hodge operator for a given cell using the COST algo. Hodge op. from primal faces to dua...
Definition cs_hodge.c:2301
void cs_hodge_epfd_cost_get(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local Hodge operator for a given cell using the COST algo. Hodge op. from primal edges to dua...
Definition cs_hodge.c:2166
void cs_hodge_edfp_cost_get(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local Hodge operator for a given cell using the COST algo. Hodge op. from dual edges to prima...
Definition cs_hodge.c:2436
void cs_hodge_vcb_wbs_get(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local Hodge operator for a given cell using the WBS algo. This function is specific for verte...
Definition cs_hodge.c:1830
void cs_hodge_fb_cost_get_stiffness(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local stiffness matrix using the generic COST algo. The computed matrix is stored in cb->loc ...
Definition cs_hodge.c:690
void cs_hodge_vcb_get_stiffness(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local stiffness matrix using the generic WBS algo. WBS standing for Whitney Barycentric Subdi...
Definition cs_hodge.c:1610
void() cs_hodge_t(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local operator for a given cell (stored in cb->hdg for a discrete Hodge operator or in cb->lo...
Definition cs_hodge.h:63
void cs_hodge_vb_cost_get_aniso_stiffness(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local stiffness matrix using the generic COST algo. Case of CDO vertex-based schemes.
Definition cs_hodge.c:989
void cs_hodge_vpcd_wbs_get(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local Hodge operator for a given cell using WBS algo. Hodge op. from primal vertices to dual ...
Definition cs_hodge.c:1946
void cs_hodge_vb_voro_get_stiffness(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local stiffness matrix using the Voronoi algorithm The computed matrix is stored in cb->loc C...
Definition cs_hodge.c:1384
void cs_hodge_vpcd_voro_get(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local Hodge operator for a given cell using VORONOI algo. Hodge op. from primal vertices to d...
Definition cs_hodge.c:2056
void cs_hodge_vb_wbs_get_stiffness(const cs_param_hodge_t h_info, const cs_cell_mesh_t *cm, cs_cell_builder_t *cb)
Build a local stiffness matrix using the generic WBS algo. WBS standing for Whitney Barycentric Subdi...
Definition cs_hodge.c:1476
Definition cs_cdo_connect.h:74
Definition cs_cdo_quantities.h:94
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
Set of local quantities and connectivities related to a mesh face Structure used to get a better memo...
Definition cs_cdo_local.h:212
Definition cs_param_cdo.h:128
Definition cs_property.h:104