File: cs_equation_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 (864 lines) | stat: -rw-r--r-- 29,814 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
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
/*============================================================================
 * Routines to handle specific settings related to a cs_equation_t structure
 *============================================================================*/

/*
  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 <assert.h>

#if defined(HAVE_PETSC)
#include <petscversion.h>
#include <petscdraw.h>
#include <petscviewer.h>
#include <petscksp.h>
#endif

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

#include <bft_error.h>
#include <bft_printf.h>
#include <bft_mem.h>

#include "cs_mesh_location.h"
#include "cs_multigrid.h"

#if defined(HAVE_PETSC)
#include "cs_sles_petsc.h"
#endif

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

#include "cs_equation_param.h"

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

BEGIN_C_DECLS

/*============================================================================
 * Type definitions
 *============================================================================*/

/*============================================================================
 * Local private variables
 *============================================================================*/

/* Default initialization */
static cs_equation_algo_t _algo_info_by_default = {
#if defined(HAVE_PETSC)
  CS_EQUATION_ALGO_PETSC_ITSOL, // Family of iterative solvers
#else
  CS_EQUATION_ALGO_CS_ITSOL,    // Family of iterative solvers
#endif
  0,                            // n_iters
  50,                           // max. number of iterations
  0,                            // n_cumulated_iters
  10000,                        // max. number of cumulated iterations
  1e-6                          // stopping criterion
};

static cs_param_itsol_t _itsol_info_by_default = {
#if defined(HAVE_PETSC)
  CS_PARAM_PRECOND_ILU0,  // preconditioner
  CS_PARAM_ITSOL_BICG,    // iterative solver
#else
  CS_PARAM_PRECOND_DIAG,  // preconditioner
  CS_PARAM_ITSOL_CG,      // iterative solver
#endif
  2500,                   // max. number of iterations
  1e-12,                  // stopping criterion on the accuracy
  150,                    // output frequency
  false                   // normalization of the residual (true or false)
};

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

#if defined(HAVE_PETSC)

/*----------------------------------------------------------------------------
 * \brief Add visualization of the matrix graph
 *
 * \param[in]  ksp     Krylov SubSpace structure
 *----------------------------------------------------------------------------*/

static inline void
_add_view(KSP          ksp)
{
  const char *p = getenv("CS_USER_PETSC_MAT_VIEW");

  if (p != NULL) {

    /* Get system and preconditioner matrixes */

    Mat a, pa;
    KSPGetOperators(ksp, &a, &pa);

    /* Output matrix in several ways depending on
       CS_USER_PETSC_MAT_VIEW environment variable */

    if (strcmp(p, "DEFAULT") == 0)
      MatView(a, PETSC_VIEWER_DEFAULT);

    else if (strcmp(p, "DRAW_WORLD") == 0)
      MatView(a, PETSC_VIEWER_DRAW_WORLD);

    else if (strcmp(p, "DRAW") == 0) {

      PetscViewer viewer;
      PetscDraw draw;
      PetscViewerDrawOpen(PETSC_COMM_WORLD, NULL, "PETSc View",
                          0, 0, 600, 600, &viewer);
      PetscViewerDrawGetDraw(viewer, 0, &draw);
      PetscViewerDrawSetPause(viewer, -1);
      MatView(a, viewer);
      PetscDrawPause(draw);

      PetscViewerDestroy(&viewer);

    }

  }

}

/*----------------------------------------------------------------------------
 * \brief Set PETSc solver and preconditioner
 *
 * \param[in, out] context  pointer to optional (untyped) value or structure
 * \param[in, out] ksp      pointer to PETSc KSP context
 *----------------------------------------------------------------------------*/

static void
_petsc_setup_hook(void   *context,
                  KSP     ksp)
{
  PC pc;

  cs_equation_param_t  *eqp = (cs_equation_param_t *)context;
  cs_param_itsol_t  info = eqp->itsol_info;
  cs_param_precond_type_t  precond = info.precond;

  /* Set the solver */
  switch (info.solver) {

  case CS_PARAM_ITSOL_CG:  /* Preconditioned Conjugate Gradient */
    KSPSetType(ksp, KSPCG);
    break;

  case CS_PARAM_ITSOL_GMRES:  /* Preconditioned GMRES */
    {
      const int  n_max_restart = 30;

      KSPSetType(ksp, KSPGMRES);
      KSPGMRESSetRestart(ksp, n_max_restart);
    }
    break;

  case CS_PARAM_ITSOL_BICG: /* Preconditioned Bi-CG */
    KSPSetType(ksp, KSPBICG);
    break;

  case CS_PARAM_ITSOL_BICGSTAB2: /* Preconditioned BiCGstab2 */
    KSPSetType(ksp, KSPBCGS);
    break;

  default:
    bft_error(__FILE__, __LINE__, 0,
              " Iterative solver not interfaced with PETSc.");
  }

  /* Set the preconditioner */
  KSPGetPC(ksp, &pc);

  if (cs_glob_n_ranks > 1) {
    if (precond == CS_PARAM_PRECOND_SSOR ||
        precond == CS_PARAM_PRECOND_ILU0) {
      precond = CS_PARAM_PRECOND_BJACOB;
      cs_base_warn(__FILE__, __LINE__);
      bft_printf(" Modify the requested preconditioner to able a parallel"
                 " computation with PETSC.\n"
                 " Switch to a block jacobi preconditioner.\n"
                 " Please check your settings.");
    }
  }


  switch (precond) {

  case CS_PARAM_PRECOND_DIAG:
    PCSetType(pc, PCJACOBI);  /* Jacobi (diagonal) preconditioning */
    break;
  case CS_PARAM_PRECOND_BJACOB:
    PCSetType(pc, PCBJACOBI);  /* Block-Jacobi (diagonal) preconditioning */
#if PETSC_VERSION_GE(3,7,0)
    PetscOptionsSetValue(NULL, "-sub_pc_factor_levels", "2");
#else
    PetscOptionsSetValue("-sub_pc_factor_levels", "2");
#endif
    break;
  case CS_PARAM_PRECOND_SSOR:
    PCSetType(pc, PCSOR);
    PCSORSetSymmetric(pc, SOR_SYMMETRIC_SWEEP);
    break;
  case CS_PARAM_PRECOND_ICC0:
    PCSetType(pc, PCICC);
    PCFactorSetLevels(pc, 0);
    break;
  case CS_PARAM_PRECOND_ILU0:
    PCSetType(pc, PCILU);
    PCFactorSetLevels(pc, 0);
    break;

  case CS_PARAM_PRECOND_AS:
    PCSetType(pc, PCASM);
    break;

  case CS_PARAM_PRECOND_AMG:
    {
      int  amg_type = 1;

      if (amg_type == 0) { // GAMG
#if PETSC_VERSION_GE(3,7,0)
        PetscOptionsSetValue(NULL, "-pc_gamg_agg_nsmooths", "1");
        PetscOptionsSetValue(NULL, "-mg_levels_ksp_type", "richardson");
        PetscOptionsSetValue(NULL, "-mg_levels_pc_type", "sor");
        PetscOptionsSetValue(NULL, "-mg_levels_ksp_max_it", "1");
        PetscOptionsSetValue(NULL, "-pc_gamg_threshold", "0.02");
        PetscOptionsSetValue(NULL, "-pc_gamg_reuse_interpolation", "TRUE");
        PetscOptionsSetValue(NULL, "-pc_gamg_square_graph", "4");
#else
        PetscOptionsSetValue("-pc_gamg_agg_nsmooths", "1");
        PetscOptionsSetValue("-mg_levels_ksp_type", "richardson");
        PetscOptionsSetValue("-mg_levels_pc_type", "sor");
        PetscOptionsSetValue("-mg_levels_ksp_max_it", "1");
        PetscOptionsSetValue("-pc_gamg_threshold", "0.02");
        PetscOptionsSetValue("-pc_gamg_reuse_interpolation", "TRUE");
        PetscOptionsSetValue("-pc_gamg_square_graph", "4");
#endif
        PCSetType(pc, PCGAMG);
      }
      else if (amg_type == 1) { // Boomer AMG (hypre)
#if PETSC_VERSION_GE(3,7,0)
        PetscOptionsSetValue(NULL,"-pc_type", "hypre");
        PetscOptionsSetValue(NULL,"-pc_hypre_type","boomeramg");
        PetscOptionsSetValue(NULL,"-pc_hypre_boomeramg_coarsen_type","HMIS");
        PetscOptionsSetValue(NULL,"-pc_hypre_boomeramg_interp_type","ext+i-cc");
        PetscOptionsSetValue(NULL,"-pc_hypre_boomeramg_agg_nl","2");
        PetscOptionsSetValue(NULL,"-pc_hypre_boomeramg_P_max","4");
        PetscOptionsSetValue(NULL,"-pc_hypre_boomeramg_strong_threshold","0.5");
        PetscOptionsSetValue(NULL,"-pc_hypre_boomeramg_no_CF","");
#else
        PetscOptionsSetValue("-pc_type", "hypre");
        PetscOptionsSetValue("-pc_hypre_type","boomeramg");
        PetscOptionsSetValue("-pc_hypre_boomeramg_coarsen_type","HMIS");
        PetscOptionsSetValue("-pc_hypre_boomeramg_interp_type","ext+i-cc");
        PetscOptionsSetValue("-pc_hypre_boomeramg_agg_nl","2");
        PetscOptionsSetValue("-pc_hypre_boomeramg_P_max","4");
        PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold","0.5");
        PetscOptionsSetValue("-pc_hypre_boomeramg_no_CF","");
#endif

        PCSetType(pc, PCHYPRE);
      }

    }
    break;

  default:
    bft_error(__FILE__, __LINE__, 0,
              " Preconditioner not interfaced with PETSc.");
  }

  /* Try to have "true" norm */
  KSPSetNormType(ksp, KSP_NORM_UNPRECONDITIONED);

  _add_view(ksp);
}

#endif /* defined(HAVE_PETSC) */

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

/*----------------------------------------------------------------------------*/
/*!
 * \brief  Create a cs_equation_param_t
 *
 * \param[in] type             type of equation
 * \param[in] var_type         type of variable (scalar, vector, tensor...)
 * \param[in] default_bc       type of boundary condition set by default
 *
 * \return a pointer to a new allocated cs_equation_param_t structure
 */
/*----------------------------------------------------------------------------*/

cs_equation_param_t *
cs_equation_param_create(cs_equation_type_t     type,
                         cs_param_var_type_t    var_type,
                         cs_param_bc_type_t     default_bc)
{
  cs_equation_param_t  *eqp = NULL;

  BFT_MALLOC(eqp, 1, cs_equation_param_t);

  eqp->type = type;
  eqp->var_type = var_type;
  eqp->verbosity = 0;
  eqp->sles_verbosity = 0;
  eqp->process_flag = 0;

  /* Build the equation flag */
  eqp->flag = 0;
  eqp->space_scheme = CS_SPACE_SCHEME_CDOVB;

  /* Vertex-based schemes imply the two following discrete Hodge operators
     Default initialization is made in accordance with this choice */
  eqp->time_hodge.inv_pty = false; // inverse property ?
  eqp->time_hodge.type = CS_PARAM_HODGE_TYPE_VPCD;
  eqp->time_hodge.algo = CS_PARAM_HODGE_ALGO_VORONOI;
  eqp->time_property = NULL;

  /* Description of the time discretization (default values) */
  eqp->time_info.scheme = CS_TIME_SCHEME_IMPLICIT;
  eqp->time_info.theta = 1.0;
  eqp->time_info.do_lumping = false;

  /* Initial condition (zero value by default) */
  eqp->time_info.n_ic_definitions = 0;
  eqp->time_info.ic_definitions = NULL;

  /* Diffusion term */
  eqp->diffusion_property = NULL;
  eqp->diffusion_hodge.inv_pty = false; // inverse property ?
  eqp->diffusion_hodge.type = CS_PARAM_HODGE_TYPE_EPFD;
  eqp->diffusion_hodge.algo = CS_PARAM_HODGE_ALGO_COST;
  eqp->diffusion_hodge.coef = 1./3.; // DGA algo.

  /* Advection term */
  eqp->advection_info.formulation = CS_PARAM_ADVECTION_FORM_CONSERV;
  eqp->advection_info.scheme = CS_PARAM_ADVECTION_SCHEME_UPWIND;
  eqp->advection_info.weight_criterion = CS_PARAM_ADVECTION_WEIGHT_XEXC;
  eqp->advection_info.quad_type = CS_QUADRATURE_BARY;
  eqp->advection_field = NULL;

  /* No reaction term by default */
  eqp->reaction_hodge.inv_pty = false;
  eqp->reaction_hodge.algo = CS_PARAM_HODGE_ALGO_WBS;
  eqp->reaction_hodge.type = CS_PARAM_HODGE_TYPE_VPCD;
  eqp->n_reaction_terms = 0;
  eqp->reaction_info = NULL;
  eqp->reaction_properties = NULL;

  /* No source term by default (always in the right-hand side) */
  eqp->n_source_terms = 0;
  eqp->source_terms = NULL;

  /* Boundary conditions structure.
     One assigns a boundary condition by default */
  eqp->bc = cs_param_bc_create(default_bc);

  /* Settings for driving the linear algebra */
  eqp->algo_info = _algo_info_by_default;
  eqp->itsol_info = _itsol_info_by_default;

  return eqp;
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief  Free a cs_equation_param_t
 *
 * \param[in, out] eqp          pointer to a cs_equation_param_t
 *
 * \return a NULL pointer
 */
/*----------------------------------------------------------------------------*/

cs_equation_param_t *
cs_equation_param_free(cs_equation_param_t     *eqp)
{
  if (eqp == NULL)
    return NULL;

  if (eqp->bc != NULL) { // Boundary conditions
    if (eqp->bc->n_defs > 0)
      BFT_FREE(eqp->bc->defs);
    BFT_FREE(eqp->bc);
    eqp->bc = NULL;
  }

  if (eqp->n_reaction_terms > 0) { // reaction terms

    for (int i = 0; i< eqp->n_reaction_terms; i++)
      BFT_FREE(eqp->reaction_info[i].name);

    BFT_FREE(eqp->reaction_info);
    BFT_FREE(eqp->reaction_properties);

    /* Remark: properties are freed when the global cs_domain_t structure is
       freed thanks to a call to cs_property_free() */
  }

  if (eqp->n_source_terms > 0) { // Source terms

    for (int i = 0; i< eqp->n_source_terms; i++)
      eqp->source_terms[i] = cs_source_term_free(eqp->source_terms[i]);
    BFT_FREE(eqp->source_terms);

  }

  cs_param_time_t  t_info = eqp->time_info;
  if (t_info.n_ic_definitions > 0) {
    for (int i = 0; i < t_info.n_ic_definitions; i++) {
      cs_param_def_t  *ic = t_info.ic_definitions + i;
      BFT_FREE(ic->ml_name);
    }
    BFT_FREE(t_info.ic_definitions);
  }

  BFT_FREE(eqp);

  return NULL;
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief  Summary of a cs_equation_param_t structure
 *
 * \param[in]  eqname   name of the related equation
 * \param[in]  eqp      pointer to a cs_equation_param_t structure
 */
/*----------------------------------------------------------------------------*/

void
cs_equation_param_summary(const char                  *eqname,
                          const cs_equation_param_t   *eqp)
{
  if (eqp == NULL)
    return;

switch (eqp->type) {
  case CS_EQUATION_TYPE_USER:
    bft_printf("  <%s/type> User-defined\n", eqname);
    break;
  case CS_EQUATION_TYPE_PREDEFINED:
    bft_printf("  <%s/type> Predefined\n", eqname);
    break;
  case CS_EQUATION_TYPE_GROUNDWATER:
    bft_printf("  <%s/type> Associated to groundwater flows\n", eqname);
    break;
  default:
    bft_error(__FILE__, __LINE__, 0,
              " Eq. %s has no type.\n"
              " Please check your settings.", eqname);
  }

  if (eqp->space_scheme == CS_SPACE_SCHEME_CDOVB)
    bft_printf("  <%s/space scheme>  CDO vertex-based\n", eqname);
  else if (eqp->space_scheme == CS_SPACE_SCHEME_CDOVCB)
    bft_printf("  <%s/space scheme>  CDO vertex+cell-based\n", eqname);
  else if (eqp->space_scheme == CS_SPACE_SCHEME_CDOFB)
    bft_printf("  <%s/space scheme>  CDO face-based\n", eqname);

  bool  unsteady = (eqp->flag & CS_EQUATION_UNSTEADY) ? true : false;
  bool  convection = (eqp->flag & CS_EQUATION_CONVECTION) ? true : false;
  bool  diffusion = (eqp->flag & CS_EQUATION_DIFFUSION) ? true : false;
  bool  reaction = (eqp->flag & CS_EQUATION_REACTION) ? true : false;
  bool  source_term = (eqp->n_source_terms > 0) ? true : false;

  bft_printf("  <%s/Terms>  unsteady:%s, convection:%s, diffusion:%s,"
             " reaction:%s, source term:%s\n",
             eqname, cs_base_strtf(unsteady), cs_base_strtf(convection),
             cs_base_strtf(diffusion), cs_base_strtf(reaction),
             cs_base_strtf(source_term));

  /* Boundary conditions */
  if (eqp->verbosity > 0) {
    cs_param_bc_t  *bcp = eqp->bc;

    bft_printf("  <%s/Boundary Conditions>\n", eqname);
    bft_printf("\t<BC/Default> %s\n", cs_param_get_bc_name(bcp->default_bc));
    if (eqp->verbosity > 1)
      bft_printf("\t<BC/Enforcement> %s\n",
                 cs_param_get_bc_enforcement_name(bcp->enforcement));
    bft_printf("\t<BC/N_Definitions> %d\n", bcp->n_defs);
    if (eqp->verbosity > 1) {
      for (int id = 0; id < bcp->n_defs; id++)
        bft_printf("\t\t<BC> Location: %s; Type: %s; Definition type: %s\n",
                   cs_mesh_location_get_name(bcp->defs[id].loc_id),
                   cs_param_get_bc_name(bcp->defs[id].bc_type),
                   cs_param_get_def_type_name(bcp->defs[id].def_type));
    }
  }

  if (unsteady) {

    const cs_param_time_t  t_info = eqp->time_info;
    const cs_param_hodge_t  h_info = eqp->time_hodge;

    bft_printf("\n  <%s/Unsteady term>\n", eqname);
    bft_printf("  <Time/Initial condition> number of definitions %d\n",
               t_info.n_ic_definitions);
    for (int i = 0; i < t_info.n_ic_definitions; i++) {
      const cs_param_def_t  *ic = t_info.ic_definitions + i;
      bft_printf("\t<Time/Initial condition> Location %s; Definition %s\n",
                 ic->ml_name, cs_param_get_def_type_name(ic->def_type));
    }
    bft_printf("  <Time/Scheme> ");
    switch (t_info.scheme) {
    case CS_TIME_SCHEME_IMPLICIT:
      bft_printf("implicit\n");
      break;
    case CS_TIME_SCHEME_EXPLICIT:
      bft_printf("explicit\n");
      break;
    case CS_TIME_SCHEME_CRANKNICO:
      bft_printf("Crank-Nicolson\n");
      break;
    case CS_TIME_SCHEME_THETA:
      bft_printf("theta scheme with value %f\n", t_info.theta);
      break;
    default:
      bft_error(__FILE__, __LINE__, 0, " Invalid time scheme.");
      break;
    }
    bft_printf("  <Time/Mass lumping> %s\n", cs_base_strtf(t_info.do_lumping));
    bft_printf("  <Time/Property> %s\n",
               cs_property_get_name(eqp->time_property));

    if (eqp->verbosity > 0) {
      bft_printf("  <Time/Hodge> %s - %s\n",
                 cs_param_hodge_get_type_name(h_info),
                 cs_param_hodge_get_algo_name(h_info));
      bft_printf("\t<Time/Hodge> Inversion of property: %s\n",
                 cs_base_strtf(h_info.inv_pty));
      if (h_info.algo == CS_PARAM_HODGE_ALGO_COST)
        bft_printf("\t<Time/Hodge> Value of the coercivity coef.: %.3e\n",
                   h_info.coef);
    }

  } /* Unsteady term */

  if (diffusion) {

    const cs_param_hodge_t  h_info = eqp->diffusion_hodge;

    bft_printf("\n  <%s/Diffusion term>\n", eqname);
    bft_printf("  <Diffusion> Property: %s\n",
               cs_property_get_name(eqp->diffusion_property));

    if (eqp->verbosity > 0) {
      bft_printf("  <Diffusion/Hodge> %s - %s\n",
                 cs_param_hodge_get_type_name(h_info),
                 cs_param_hodge_get_algo_name(h_info));
      bft_printf("\t<Diffusion/Hodge> Inversion of property: %s\n",
                 cs_base_strtf(h_info.inv_pty));
      if (h_info.algo == CS_PARAM_HODGE_ALGO_COST)
        bft_printf("\t<Diffusion/Hodge> Value of the coercivity coef.: %.3e\n",
                   h_info.coef);
    }

  } /* Diffusion term */

  if (convection) {

    const cs_param_advection_t  a_info = eqp->advection_info;

    bft_printf("\n  <%s/Advection term>\n", eqname);
    bft_printf("  <Advection field>  %s\n",
               cs_advection_field_get_name(eqp->advection_field));

    if (eqp->verbosity > 0) {
      bft_printf("  <Advection/Formulation>");
      switch(a_info.formulation) {
      case CS_PARAM_ADVECTION_FORM_CONSERV:
        bft_printf(" Conservative\n");
        break;
      case CS_PARAM_ADVECTION_FORM_NONCONS:
        bft_printf(" Non-conservative\n");
        break;
      default:
        bft_error(__FILE__, __LINE__, 0,
                  " Invalid operator type for advection.");
      }

      bft_printf("  <Advection/Scheme> ");
      switch(a_info.scheme) {
      case CS_PARAM_ADVECTION_SCHEME_CENTERED:
        bft_printf(" centered\n");
        break;
      case CS_PARAM_ADVECTION_SCHEME_CIP:
        bft_printf(" continuous interior penalty\n");
        break;
      case CS_PARAM_ADVECTION_SCHEME_UPWIND:
        bft_printf(" upwind\n");
        break;
      case CS_PARAM_ADVECTION_SCHEME_SAMARSKII:
        bft_printf(" upwind weighted with Samarskii function\n");
        break;
      case CS_PARAM_ADVECTION_SCHEME_SG:
        bft_printf(" upwind weighted with Scharfetter-Gummel function\n");
        break;
      default:
        bft_error(__FILE__, __LINE__, 0,
                  " Invalid weight algorithm for advection.");
      }

    } // verbosity > 0

  } // Advection term

  if (reaction) {

    bft_printf("\n  <%s/Number of reaction terms> %d\n",
               eqname, eqp->n_reaction_terms);

    if (eqp->verbosity > 0) {

      const cs_param_hodge_t  h_info = eqp->reaction_hodge;

      bft_printf("  <Reaction/Hodge> %s - %s\n",
                 cs_param_hodge_get_type_name(h_info),
                 cs_param_hodge_get_algo_name(h_info));
      if (h_info.algo == CS_PARAM_HODGE_ALGO_COST)
        bft_printf("\t<Reaction/Hodge> Value of the coercivity coef.: %.3e\n",
                   h_info.coef);

    }

    for (int r_id = 0; r_id < eqp->n_reaction_terms; r_id++) {
      cs_param_reaction_t  r_info = eqp->reaction_info[r_id];
      bft_printf("  <Reaction tem %02d> Property: %s; Operator type: %s\n",
                 r_id, cs_property_get_name(eqp->reaction_properties[r_id]),
                 cs_param_reaction_get_type_name(r_info.type));

    }

  } // Reaction terms

  if (source_term) {

    bft_printf("\n  <%s/Source terms>\n", eqname);
    for (int s_id = 0; s_id < eqp->n_source_terms; s_id++)
      cs_source_term_summary(eqname, eqp->source_terms[s_id]);

  } // Source terms

  /* Iterative solver information */
  const cs_param_itsol_t   itsol = eqp->itsol_info;

  bft_printf("\n  <%s/Sparse Linear Algebra>", eqname);
  if (eqp->algo_info.type == CS_EQUATION_ALGO_CS_ITSOL)
    bft_printf(" Code_Saturne iterative solvers\n");
  else if (eqp->algo_info.type == CS_EQUATION_ALGO_PETSC_ITSOL)
    bft_printf(" PETSc iterative solvers\n");
  bft_printf("\t<sla> Solver.MaxIter     %d\n", itsol.n_max_iter);
  bft_printf("\t<sla> Solver.Name        %s\n",
             cs_param_get_solver_name(itsol.solver));
  bft_printf("\t<sla> Solver.Precond     %s\n",
             cs_param_get_precond_name(itsol.precond));
  bft_printf("\t<sla> Solver.Eps        % -10.6e\n", itsol.eps);
  bft_printf("\t<sla> Solver.Normalized  %s\n",
             cs_base_strtf(itsol.resid_normalized));

}

/*----------------------------------------------------------------------------*/
/*!
 * \brief Initialize SLES structure for the resolution of the linear system
 *        according to the settings related to this equation
 *
 * \param[in]   eqname       pointer to an cs_equation_t structure
 * \param[in]   eqp          pointer to a cs_equation_param_t struct.
 * \param[in]   field_id     id of the cs_field_t struct. for this equation
 */
/*----------------------------------------------------------------------------*/

void
cs_equation_param_init_sles(const char                 *eqname,
                            const cs_equation_param_t  *eqp,
                            int                         field_id)
{

  const cs_equation_algo_t  algo = eqp->algo_info;
  const cs_param_itsol_t  itsol = eqp->itsol_info;

  switch (algo.type) {
  case CS_EQUATION_ALGO_CS_ITSOL:
    {
      int  poly_degree = 0; // by default: Jacobi preconditioner

      if (itsol.precond == CS_PARAM_PRECOND_POLY1)
        poly_degree = 1;

      if (itsol.precond != CS_PARAM_PRECOND_POLY1 &&
          itsol.precond != CS_PARAM_PRECOND_DIAG)
        bft_error(__FILE__, __LINE__, 0,
                  " Incompatible preconditioner with Code_Saturne solvers.\n"
                  " Please change your settings (try PETSc ?)");

      switch (itsol.solver) { // Type of iterative solver
      case CS_PARAM_ITSOL_CG:
        cs_sles_it_define(field_id,  // give the field id (future: eq_id ?)
                          NULL,
                          CS_SLES_PCG,
                          poly_degree,
                          itsol.n_max_iter);
        break;
      case CS_PARAM_ITSOL_BICG:
        cs_sles_it_define(field_id,  // give the field id (future: eq_id ?)
                          NULL,
                          CS_SLES_BICGSTAB,
                          poly_degree,
                          itsol.n_max_iter);
        break;
      case CS_PARAM_ITSOL_BICGSTAB2:
        cs_sles_it_define(field_id,  // give the field id (future: eq_id ?)
                          NULL,
                          CS_SLES_BICGSTAB2,
                          poly_degree,
                          itsol.n_max_iter);
        break;
      case CS_PARAM_ITSOL_CR3:
        cs_sles_it_define(field_id,  // give the field id (future: eq_id ?)
                          NULL,
                          CS_SLES_PCR3,
                          poly_degree,
                          itsol.n_max_iter);
        break;
      case CS_PARAM_ITSOL_GMRES:
        cs_sles_it_define(field_id,  // give the field id (future: eq_id ?)
                          NULL,
                          CS_SLES_GMRES,
                          poly_degree,
                          itsol.n_max_iter);
        break;
      case CS_PARAM_ITSOL_AMG:
        {
          cs_multigrid_t  *mg = cs_multigrid_define(field_id, NULL);

          /* Advanced setup (default is specified inside the brackets) */
          cs_multigrid_set_solver_options
            (mg,
             CS_SLES_JACOBI,   // descent smoother type (CS_SLES_PCG)
             CS_SLES_JACOBI,   // ascent smoother type (CS_SLES_PCG)
             CS_SLES_PCG,      // coarse solver type (CS_SLES_PCG)
             itsol.n_max_iter, // n max cycles (100)
             5,                // n max iter for descent (10)
             5,                // n max iter for asscent (10)
             1000,             // n max iter coarse solver (10000)
             0,                // polynomial precond. degree descent (0)
             0,                // polynomial precond. degree ascent (0)
             0,                // polynomial precond. degree coarse (0)
             1.0,    // precision multiplier descent (< 0 forces max iters)
             1.0,    // precision multiplier ascent (< 0 forces max iters)
             1);     // requested precision multiplier coarse (default 1)

        }
      default:
        bft_error(__FILE__, __LINE__, 0,
                  _(" Undefined iterative solver for solving %s equation.\n"
                    " Please modify your settings."), eqname);
        break;
      } // end of switch

      /* Define the level of verbosity for SLES structure */
      if (eqp->sles_verbosity > 1) {

        cs_sles_t  *sles = cs_sles_find_or_add(field_id, NULL);
        cs_sles_it_t  *sles_it = (cs_sles_it_t *)cs_sles_get_context(sles);

        /* Set verbosity */
        cs_sles_set_verbosity(sles, eqp->sles_verbosity);

        if (eqp->sles_verbosity > 2) /* Add plot */
          cs_sles_it_set_plot_options(sles_it, eqname,
                                      true);    /* use_iteration instead of
                                                   wall clock time */

      }

    } // Solver provided by Code_Saturne
    break;

  case CS_EQUATION_ALGO_PETSC_ITSOL:
    {
#if defined(HAVE_PETSC)

      /* Initialization must be called before setting options;
         it does not need to be called before calling
         cs_sles_petsc_define(), as this is handled automatically. */

      PetscBool is_initialized;
      PetscInitialized(&is_initialized);
      if (is_initialized == PETSC_FALSE) {
#if defined(HAVE_MPI)
        PETSC_COMM_WORLD = cs_glob_mpi_comm;
#endif
        PetscInitializeNoArguments();
      }

      if (eqp->itsol_info.precond == CS_PARAM_PRECOND_SSOR ||
          eqp->itsol_info.precond == CS_PARAM_PRECOND_ILU0 ||
          eqp->itsol_info.precond == CS_PARAM_PRECOND_ICC0)
        cs_sles_petsc_define(field_id,
                             NULL,
                             MATSEQAIJ, // Warning SEQ not MPI
                             _petsc_setup_hook,
                             (void *)eqp);
      else
        cs_sles_petsc_define(field_id,
                             NULL,
                             MATMPIAIJ,
                             _petsc_setup_hook,
                             (void *)eqp);
#else
      bft_error(__FILE__, __LINE__, 0,
                _(" PETSC algorithms used to solve %s are not linked.\n"
                  " Please install Code_Saturne with PETSc."), eqname);

#endif // HAVE_PETSC
    } // Solver provided by PETSc
    break;

  default:
    bft_error(__FILE__, __LINE__, 0,
              _(" Algorithm requested to solve %s is not implemented yet.\n"
                " Please modify your settings."), eqname);
    break;

  } // end switch on algorithms

}