File: cs_param.c

package info (click to toggle)
code-saturne 4.3.3%2Brepack-1
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 77,992 kB
  • sloc: ansic: 281,257; f90: 122,305; python: 56,490; makefile: 3,915; xml: 3,285; cpp: 3,183; sh: 1,139; lex: 176; yacc: 101; sed: 16
file content (604 lines) | stat: -rw-r--r-- 17,680 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
/*============================================================================
 * Routines to handle the definition and usage of material properties
 *============================================================================*/

/*
  This file is part of Code_Saturne, a general-purpose CFD tool.

  Copyright (C) 1998-2016 EDF S.A.

  This program is free software; you can redistribute it and/or modify it under
  the terms of the GNU General Public License as published by the Free Software
  Foundation; either version 2 of the License, or (at your option) any later
  version.

  This program is distributed in the hope that it will be useful, but WITHOUT
  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
  FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
  details.

  You should have received a copy of the GNU General Public License along with
  this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
  Street, Fifth Floor, Boston, MA 02110-1301, USA.
*/

/*----------------------------------------------------------------------------*/

#include "cs_defs.h"

/*----------------------------------------------------------------------------
 * Standard C library headers
 *----------------------------------------------------------------------------*/

#include <stdlib.h>
#include <stdio.h>
#include <string.h>
#include <assert.h>

/*----------------------------------------------------------------------------
 *  Local headers
 *----------------------------------------------------------------------------*/

#include "bft_mem.h"
#include "bft_printf.h"

#include "cs_mesh_location.h"
#include "cs_field.h"
#include "cs_cdo.h"

/*----------------------------------------------------------------------------
 * Header for the current file
 *----------------------------------------------------------------------------*/

#include "cs_param.h"

/*----------------------------------------------------------------------------*/

BEGIN_C_DECLS

/*=============================================================================
 * Local Macro definitions and structure definitions
 *============================================================================*/

#define CS_PARAM_DBG  0

/*============================================================================
 * Global static variables
 *============================================================================*/

static const char
cs_param_def_type_name[CS_PARAM_N_DEF_TYPES][CS_BASE_STRING_LEN]=
  { N_("by analytic function"),
    N_("by array"),
    N_("by law (one argument)"),
    N_("by law (two arguments)"),
    N_("by law (two arguments: scalar+vector)"),
    N_("quantity over a volume"),
    N_("by time function"),
    N_("by user function"),
    N_("by value") };

static const char
cs_param_var_type_name[CS_PARAM_N_VAR_TYPES][CS_BASE_STRING_LEN]=
  { N_("scalar"),
    N_("vector"),
    N_("tensor") };

static const char
cs_param_bc_type_name[CS_PARAM_N_BC_TYPES][CS_BASE_STRING_LEN] =
  { N_("Homogeneous Dirichlet"),
    N_("Dirichlet"),
    N_("Homogeneous Neumann"),
    N_("Neumann"),
    N_("Robin") };

static const char
cs_param_boundary_type_name[CS_PARAM_N_BOUNDARY_TYPES][CS_BASE_STRING_LEN] =
  { N_("Wall"),
    N_("Inlet"),
    N_("Outlet"),
    N_("Symmetry") };

static const char
cs_param_hodge_type_desc[CS_PARAM_N_HODGE_TYPES][CS_BASE_STRING_LEN] =
  { "VpCd",
    "EpFd",
    "FpEd",
    "EdFp",
    "CpVd"  };

static const char
cs_param_hodge_algo_desc[CS_PARAM_N_HODGE_ALGOS][CS_BASE_STRING_LEN] =
  { "Voronoi",
    "Whitney on the Barycentric Subdivision (WBS)",
    "COnsistency-STabilization splitting (COST)" };

/*============================================================================
 * Private function prototypes
 *============================================================================*/

/*============================================================================
 * Public function prototypes
 *============================================================================*/

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Get the name related to a type of variable
 *
 * \param[in] type     cs_param_var_type_t
 *
 * \return the name associated to this type
 */
/*----------------------------------------------------------------------------*/

const char *
cs_param_get_var_type_name(const cs_param_var_type_t   type)
{
  return cs_param_var_type_name[type];
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Get the name related to a type of definition
 *
 * \param[in] type     cs_param_def_type_t
 *
 * \return the name associated to this type
 */
/*----------------------------------------------------------------------------*/

const char *
cs_param_get_def_type_name(const cs_param_def_type_t   type)
{
  return cs_param_def_type_name[type];
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief  Set a cs_def_t structure
 *
 * \param[in]      def_type   type of definition (by value, function...)
 * \param[in]      var_type   type of variables (scalar, vector, tensor...)
 * \param[in]      val        value to set
 * \param[in, out] def        pointer to cs_def_t structure
 */
/*----------------------------------------------------------------------------*/

void
cs_param_set_def(cs_param_def_type_t      def_type,
                 cs_param_var_type_t      var_type,
                 const void              *val,
                 cs_def_t                *def)
{
  assert(var_type != CS_PARAM_N_VAR_TYPES);

  switch (def_type) {

  case CS_PARAM_DEF_BY_VALUE:
  case CS_PARAM_DEF_BY_QOV:
    if (val == NULL) {
      if (var_type == CS_PARAM_VAR_SCAL)
        def->get.val = 0.0;
      else if (var_type == CS_PARAM_VAR_VECT)
        def->get.vect[0] = def->get.vect[1] = def->get.vect[2] = 0.0;
      else if (var_type == CS_PARAM_VAR_TENS) {
        for (int i = 0; i < 3; i++)
          for (int j = 0; j < 3; j++)
            def->get.tens[i][j] = 0.0;
      }
      else {
        assert(var_type == CS_PARAM_VAR_SYMTENS);
        for (int i = 0; i < 6; i++)
          def->get.twovects[i] = 0.0;
      }
    }
    else { // val != NULL
      if (var_type == CS_PARAM_VAR_SCAL)
        def->get.val = atof(val);
      else if (var_type == CS_PARAM_VAR_VECT) {
        char s[3][32];
        sscanf(val, "%s %s %s", s[0], s[1], s[2]);
        for (int i = 0; i < 3; i++)
          def->get.vect[i] = atof(s[i]);
      }
      else if (var_type == CS_PARAM_VAR_TENS) {
        char s[9][32];
        sscanf(val, "%s %s %s %s %s %s %s %s %s",
               s[0], s[1], s[2], s[3], s[4], s[5], s[6], s[7], s[8]);
        for (int i = 0; i < 3; i++)
          for (int j = 0; j < 3; j++)
            def->get.tens[i][j] = atof(s[3*i+j]);
      }
      else {
        assert(var_type == CS_PARAM_VAR_SYMTENS);
        char s[6][32];
        sscanf(val, "%s %s %s %s %s %s", s[0], s[1], s[2], s[3], s[4], s[5]);
        for (int i = 0; i < 6; i++)
          def->get.twovects[i] = atof(s[i]);
      }
    }
    break;

  case CS_PARAM_DEF_BY_ANALYTIC_FUNCTION:
    if (val == NULL)
      def->analytic = NULL;
    else
      def->analytic = (cs_analytic_func_t *)val;
    break;

  case CS_PARAM_DEF_BY_TIME_FUNCTION:
    if (val == NULL)
      def->time_func = NULL;
    else
      def->time_func = (cs_timestep_func_t *)val;
    break;

  case CS_PARAM_DEF_BY_LAW_ONESCA:
    if (val == NULL)
      def->law1_func = NULL;
    else
      def->law1_func = (cs_onevar_law_func_t *)val;
    break;

  default:
    bft_error(__FILE__, __LINE__, 0,
              " This type of definition is not handled yet.\n"
              " Please modify your settings.");
    break;

  } /* end of switch on def_type */

}

/*----------------------------------------------------------------------------*/
/*!
 * \brief  Set a cs_get_t structure
 *
 * \param[in]      var_type   type of variables (scalar, vector, tensor...)
 * \param[in]      val        value to set
 * \param[in, out] get        pointer to cs_get_t structure
 */
/*----------------------------------------------------------------------------*/

void
cs_param_set_get(cs_param_var_type_t      var_type,
                 const char              *val,
                 cs_get_t                *get)
{
  assert(var_type != CS_PARAM_N_VAR_TYPES);

  if (val == NULL) {

    if (var_type == CS_PARAM_VAR_SCAL)
      get->val = 0.0;
    else if (var_type == CS_PARAM_VAR_VECT)
      get->vect[0] = get->vect[1] = get->vect[2] = 0.0;
    else if (var_type == CS_PARAM_VAR_TENS) {
      for (int i = 0; i < 3; i++)
        for (int j = 0; j < 3; j++)
          get->tens[i][j] = 0.0;
    }
    else if (var_type == CS_PARAM_VAR_SYMTENS) {
      for (int i = 0; i < 6; i++)
        get->twovects[i] = 0.0;
    }
    else
      bft_error(__FILE__, __LINE__, 0, _(" Invalid type of variable."));

  }
  else { // val != NULL

    if (var_type == CS_PARAM_VAR_SCAL)
      get->val = atof(val);
    else if (var_type == CS_PARAM_VAR_VECT) {
      char s[3][32];
      sscanf(val, "%s %s %s", s[0], s[1], s[2]);
      for (int i = 0; i < 3; i++)
        get->vect[i] = atof(s[i]);
    }
    else if (var_type == CS_PARAM_VAR_TENS) {
      char s[9][32];
      sscanf(val, "%s %s %s %s %s %s %s %s %s",
             s[0], s[1], s[2], s[3], s[4], s[5], s[6], s[7], s[8]);
      for (int i = 0; i < 3; i++)
        for (int j = 0; j < 3; j++)
          get->tens[i][j] = atof(s[3*i+j]);
    }
    else if (var_type == CS_PARAM_VAR_SYMTENS) {
      char s[6][32];
      sscanf(val, "%s %s %s %s %s %s", s[0], s[1], s[2], s[3], s[4], s[5]);
      for (int i = 0; i < 6; i++)
        get->twovects[i] = atof(s[i]);
    }
    else
      bft_error(__FILE__, __LINE__, 0, _(" Invalid type of variable."));

  }

}

/*----------------------------------------------------------------------------*/
/*!
 * \brief  Allocate and initialize a new cs_param_bc_t structure
 *
 * \param[in]  default_bc     default boundary condition
 *
 * \return a pointer to the new structure (free with cs_equation_param_t)
 */
/*----------------------------------------------------------------------------*/

cs_param_bc_t *
cs_param_bc_create(cs_param_bc_type_t  default_bc)
{
  cs_param_bc_t  *bc = NULL;

  BFT_MALLOC(bc, 1, cs_param_bc_t);

  bc->default_bc = default_bc;
  /* Initialization by default */
  bc->enforcement = CS_PARAM_BC_ENFORCE_STRONG;
  bc->quad_type = CS_QUADRATURE_BARY;
  bc->use_subdiv = false;

  bc->n_defs = 0;
  bc->defs = NULL;

  return bc;
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief  Set a cs_param_bc_def_t structure
 *
 * \param[in, out] bcpd       pointer to cs_param_bc_def_t struct. to set
 * \param[in]      loc_id     id related to a cs_mesh_location_t
 * \param[in]      bc_type    generic type of admissible boundary conditions
 * \param[in]      var_type   type of variables (scalar, vector, tensor...)
 * \param[in]      def_type   by value, function...
 * \param[in]      coef1      access to the value of the first coef
 * \param[in]      coef2      access to the value of the second coef (optional)
 */
/*----------------------------------------------------------------------------*/

void
cs_param_bc_def_set(cs_param_bc_def_t      *bcpd,
                    int                     loc_id,
                    cs_param_bc_type_t      bc_type,
                    cs_param_var_type_t     var_type,
                    cs_param_def_type_t     def_type,
                    const void             *coef1,
                    const void             *coef2)
{
  if (bcpd == NULL)
    return;

  /* Sanity checks */
  assert(def_type != CS_PARAM_N_DEF_TYPES);
  assert(bc_type != CS_PARAM_N_BC_TYPES);

  bcpd->loc_id = loc_id;
  bcpd->var_type = var_type;
  bcpd->bc_type = bc_type;
  bcpd->def_type = def_type;

  cs_param_set_def(def_type, var_type, coef1, &(bcpd->def_coef1));
  cs_param_set_def(def_type, var_type, coef2, &(bcpd->def_coef2));
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Get the name of the type of boundary condition
 *
 * \param[in] bc_type     type of boundary condition
 *
 * \return the associated bc name
 */
/*----------------------------------------------------------------------------*/

const char *
cs_param_get_bc_name(cs_param_bc_type_t  bc)
{
  switch(bc) {

  case CS_PARAM_BC_HMG_DIRICHLET:
    return "Homogeneous Dirichlet";
    break;
  case CS_PARAM_BC_DIRICHLET:
    return "Dirichlet";
    break;
  case CS_PARAM_BC_HMG_NEUMANN:
    return "Homogeneous Neumann";
  case CS_PARAM_BC_NEUMANN:
    return "Neumann";
  case CS_PARAM_BC_ROBIN:
    return "Robin";
  default:
    bft_error(__FILE__, __LINE__, 0,
              _(" Invalid BC type. Stop execution."));
  }

  return "NULL"; // avoid a warning
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Get the name of the type of enforcement of the boundary condition
 *
 * \param[in] bc_enforce    type of enforcement of boundary conditions
 *
 * \return the associated name
 */
/*----------------------------------------------------------------------------*/

const char *
cs_param_get_bc_enforcement_name(cs_param_bc_enforce_t  type)
{
  switch(type) {

  case CS_PARAM_BC_ENFORCE_STRONG:
    return "strong";
    break;
  case CS_PARAM_BC_ENFORCE_WEAK_PENA:
    return "weak with a big penalization coefficient";
    break;
  case CS_PARAM_BC_ENFORCE_WEAK_NITSCHE:
    return "weak using the Nitsche method";
  case CS_PARAM_BC_ENFORCE_WEAK_SYM:
    return "weak using the symmetrized Nitsche method";
  default:
    bft_error(__FILE__, __LINE__, 0,
              _(" Invalid type of enforcement. Stop execution."));
  }

  return "NULL"; // avoid a warning
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Get the name of the type of reaction term
 *
 * \param[in] r_type     type of reaction term
 *
 * \return the name associated with this type of reaction term
 */
/*----------------------------------------------------------------------------*/

const char *
cs_param_reaction_get_type_name(cs_param_reaction_type_t  r_type)
{
  switch (r_type) {
  case CS_PARAM_REACTION_TYPE_LINEAR:
    return "Linear";
    break;
  case CS_PARAM_N_REACTION_TYPES:
    return "Not set";
    break;
  default:
    bft_error(__FILE__, __LINE__, 0,
              _(" Invalid type of reaction term. Stop execution."));
  }

  return "NULL"; // avoid a warning
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Get the name of algorithm related to a discrete Hdoge operator
 *
 * \param[in] h_info     cs_param_hodge_t structure
 *
 * \return the name of the algorithm
 */
/*----------------------------------------------------------------------------*/

const char *
cs_param_hodge_get_algo_name(const cs_param_hodge_t   h_info)
{
  return cs_param_hodge_algo_desc[h_info.algo];
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Get the type of discrete Hodge operator
 *
 * \param[in] h_info     cs_param_hodge_t structure
 *
 * \return the name of the type
 */
/*----------------------------------------------------------------------------*/

const char *
cs_param_hodge_get_type_name(const cs_param_hodge_t   h_info)
{
  return cs_param_hodge_type_desc[h_info.type];
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Get the name of the solver
 *
 * \param[in] solver     type of iterative solver
 *
 * \return the associated solver name
 */
/*----------------------------------------------------------------------------*/

const char *
cs_param_get_solver_name(cs_param_itsol_type_t  solver)
{
  switch (solver) {
  case CS_PARAM_ITSOL_CG:
    return  "CG";
    break;
  case CS_PARAM_ITSOL_BICG:
    return "BiCG";
    break;
  case CS_PARAM_ITSOL_BICGSTAB2:
    return "BiCGstab2";
    break;
  case CS_PARAM_ITSOL_CR3:
    return "Conjugate.Residual.3Layers";
    break;
  case CS_PARAM_ITSOL_GMRES:
    return "GMRES";
    break;
  case CS_PARAM_ITSOL_AMG:
    return "Algebraic.Multigrid";
    break;
  default:
    bft_error(__FILE__, __LINE__, 0,
              _(" Invalid solver. Stop execution."));
  }

  return "NULL";
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Get the name of the preconditioner
 *
 * \param[in] precond     type of preconditioner
 *
 * \return the associated preconditioner name
 */
/*----------------------------------------------------------------------------*/

const char *
cs_param_get_precond_name(cs_param_precond_type_t  precond)
{
  switch (precond) {
  case CS_PARAM_PRECOND_DIAG:
    return  "Diagonal";
    break;
  case CS_PARAM_PRECOND_BJACOB:
    return  "Block-Jacobi";
    break;
  case CS_PARAM_PRECOND_POLY1:
    return  "Neumann.Poly.O1";
    break;
  case CS_PARAM_PRECOND_SSOR:
    return  "SSOR";
    break;
  case CS_PARAM_PRECOND_ILU0:
    return  "ILU0";
    break;
  case CS_PARAM_PRECOND_ICC0:
    return  "ICC0";
    break;
  case CS_PARAM_PRECOND_AMG:
    return  "Algebraic.MultiGrid";
    break;
  case CS_PARAM_PRECOND_AS:
    return  "Additive.Schwarz";
    break;
  default:
    bft_error(__FILE__, __LINE__, 0,
              _(" Invalid preconditioner. Stop execution."));
  }

  return "NULL";
}

/*----------------------------------------------------------------------------*/

END_C_DECLS