1 /*============================================================================
2  * Functions associated with code coupling.
3  *============================================================================*/
4 
5 /*
6   This file is part of Code_Saturne, a general-purpose CFD tool.
7 
8   Copyright (C) 1998-2021 EDF S.A.
9 
10   This program is free software; you can redistribute it and/or modify it under
11   the terms of the GNU General Public License as published by the Free Software
12   Foundation; either version 2 of the License, or (at your option) any later
13   version.
14 
15   This program is distributed in the hope that it will be useful, but WITHOUT
16   ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
17   FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
18   details.
19 
20   You should have received a copy of the GNU General Public License along with
21   this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
22   Street, Fifth Floor, Boston, MA 02110-1301, USA.
23 */
24 
25 /*----------------------------------------------------------------------------*/
26 
27 #include "cs_defs.h"
28 
29 /*----------------------------------------------------------------------------
30  * Standard C library headers
31  *----------------------------------------------------------------------------*/
32 
33 #include <assert.h>
34 #include <math.h>
35 #include <stdarg.h>
36 #include <stdio.h>
37 #include <stdlib.h>
38 #include <string.h>
39 
40 /*----------------------------------------------------------------------------
41  * PLE library headers
42  *----------------------------------------------------------------------------*/
43 
44 #include <ple_coupling.h>
45 #include <ple_locator.h>
46 
47 /*----------------------------------------------------------------------------
48  * Local headers
49  *----------------------------------------------------------------------------*/
50 
51 #include "bft_mem.h"
52 #include "bft_printf.h"
53 
54 #include "fvm_nodal.h"
55 #include "fvm_writer.h"
56 
57 #include "cs_base.h"
58 #include "cs_coupling.h"
59 #include "cs_mesh.h"
60 #include "cs_mesh_quantities.h"
61 #include "cs_mesh_connect.h"
62 #include "cs_prototypes.h"
63 #include "cs_selector.h"
64 
65 /*----------------------------------------------------------------------------
66  *  Header for the current file
67  *----------------------------------------------------------------------------*/
68 
69 #include "cs_sat_coupling.h"
70 
71 /*----------------------------------------------------------------------------*/
72 
73 BEGIN_C_DECLS
74 
75 /*! \cond DOXYGEN_SHOULD_SKIP_THIS */
76 
77 /*=============================================================================
78  * Local Structure Definitions
79  *============================================================================*/
80 
81 /* Structure associated with Code_Saturne coupling */
82 
83 typedef struct {
84 
85   int      match_id;        /* Id of matched application, -1 initially */
86   char    *app_name;        /* Application name */
87   char    *face_cpl_sel_c;  /* Face selection criteria */
88   char    *cell_cpl_sel_c;  /* Cell selection criteria */
89   char    *face_loc_sel_c;  /* Face selection criteria */
90   char    *cell_loc_sel_c;  /* Cell selection criteria */
91   int      verbosity;       /* Verbosity level */
92 
93 } _cs_sat_coupling_builder_t;
94 
95 struct _cs_sat_coupling_t {
96 
97   char                   *sat_name;        /* Application name */
98   cs_sat_coupling_tag_t  *tag_func;        /* Tagging function pointer */
99   void                   *tag_context;     /* Tagging context */
100 
101   char            *face_cpl_sel; /* Face selection criteria */
102   char            *cell_cpl_sel; /* Face selection criteria */
103   char            *face_loc_sel; /* Face selection criteria */
104   char            *cell_loc_sel; /* Face selection criteria */
105 
106   ple_locator_t   *localis_cel;  /* Locator associated with cells */
107   ple_locator_t   *localis_fbr;  /* Locator associated with boundary faces */
108 
109   cs_lnum_t        nbr_cel_sup;  /* Number of associated cell locations */
110   cs_lnum_t        nbr_fbr_sup;  /* Number of associated face locations */
111 
112   fvm_nodal_t     *cells_sup;    /* Local cells at which distant values are
113                                     interpolated*/
114   fvm_nodal_t     *faces_sup;    /* Local faces at which distant values are
115                                     interpolated*/
116 
117   cs_real_t       *distant_dist_fbr; /* Distant vectors (distance JJ') */
118   cs_real_t       *distant_of;
119   cs_real_t       *local_of;
120   cs_real_t       *distant_pond_fbr; /* Distant weighting coefficient */
121   cs_real_t       *local_pond_fbr;   /* Local weighting coefficient */
122 
123   cs_real_t        tolerance; /* location tolerance */
124   int              verbosity; /* Verbosity level */
125 
126   /* Communication-related members */
127 
128 #if defined(HAVE_MPI)
129 
130   MPI_Comm         comm;           /* Associated MPI communicator */
131 
132   int              n_sat_ranks;    /* Number of associated Code_Saturne ranks */
133   int              sat_root_rank;  /* First associated Code_Saturne rank */
134 
135 #endif
136 
137 };
138 
139 /*============================================================================
140  * Static global variables
141  *============================================================================*/
142 
143 /* Array of couplings */
144 
145 static int _cs_glob_n_sat_cp = -1;
146 
147 static int                         _sat_coupling_builder_size = 0;
148 static _cs_sat_coupling_builder_t *_sat_coupling_builder = NULL;
149 
150 static int                  cs_glob_sat_n_couplings = 0;
151 static cs_sat_coupling_t  **cs_glob_sat_couplings = NULL;
152 
153 /*============================================================================
154  * Private function definitions
155  *============================================================================*/
156 
157 /*----------------------------------------------------------------------------
158  * Remove matched builder entries from the coupling builder.
159  *----------------------------------------------------------------------------*/
160 
161 static void
_remove_matched_builder_entries(void)162 _remove_matched_builder_entries(void)
163 {
164   int i;
165   int n_unmatched_entries = 0;
166 
167   /* First, free arrays associated with marked entries */
168 
169   for (i = 0; i < _sat_coupling_builder_size; i++) {
170 
171     _cs_sat_coupling_builder_t *scb = _sat_coupling_builder + i;
172 
173     if (scb->match_id != -1) {
174       if (scb->face_cpl_sel_c != NULL) BFT_FREE(scb->face_cpl_sel_c);
175       if (scb->cell_cpl_sel_c != NULL) BFT_FREE(scb->cell_cpl_sel_c);
176       if (scb->face_loc_sel_c != NULL) BFT_FREE(scb->face_loc_sel_c);
177       if (scb->cell_loc_sel_c != NULL) BFT_FREE(scb->cell_loc_sel_c);
178       if (scb->app_name != NULL) BFT_FREE(scb->app_name);
179     }
180   }
181 
182   /* Now, remove marked entries and resize */
183 
184   for (i = 0; i < _sat_coupling_builder_size; i++) {
185     _cs_sat_coupling_builder_t *scb = _sat_coupling_builder + i;
186     if (scb->match_id < 0) {
187       *(_sat_coupling_builder + n_unmatched_entries) = *scb;
188       n_unmatched_entries += 1;
189     }
190   }
191 
192   _sat_coupling_builder_size = n_unmatched_entries;
193 
194   BFT_REALLOC(_sat_coupling_builder,
195               _sat_coupling_builder_size,
196               _cs_sat_coupling_builder_t);
197 }
198 
199 /*----------------------------------------------------------------------------
200  * Print information on yet unmatched Code_Saturne couplings.
201  *----------------------------------------------------------------------------*/
202 
203 static void
_print_all_unmatched_sat(void)204 _print_all_unmatched_sat(void)
205 {
206   int i;
207 
208   const char empty_string[] = "";
209 
210   /* Loop on defined Code_Saturne instances */
211 
212   for (i = 0; i < _sat_coupling_builder_size; i++) {
213 
214     _cs_sat_coupling_builder_t *scb = _sat_coupling_builder + i;
215 
216     if (scb->match_id < 0) {
217 
218       const char *local_name = empty_string;
219 
220       if (scb->app_name != NULL)
221         local_name = scb->app_name;
222 
223       bft_printf(_(" Code_Saturne coupling:\n"
224                    "   coupling id:              %d\n"
225                    "   local name:               \"%s\"\n\n"),
226                  i, local_name);
227     }
228   }
229 
230   bft_printf_flush();
231 }
232 
233 /*----------------------------------------------------------------------------
234  * Initialize communicator for Code_Saturne coupling
235  *
236  * parameters:
237  *   sat_coupling  <-> Code_Saturne coupling structure
238  *   coupling_id   <-- id of this coupling (for log file message)
239  *----------------------------------------------------------------------------*/
240 
241 static void
_init_comm(cs_sat_coupling_t * sat_coupling,int coupling_id)242 _init_comm(cs_sat_coupling_t *sat_coupling,
243            int                coupling_id)
244 
245 {
246 #if defined(HAVE_MPI)
247 
248   int  mpi_flag = 0;
249   int local_range[2] = {-1, -1};
250   int distant_range[2] = {-1, -1};
251 
252   MPI_Initialized(&mpi_flag);
253 
254   if (mpi_flag == 0)
255     return;
256 
257   bft_printf(_(" Code_Saturne coupling %d: initializing MPI communication ... "),
258              coupling_id);
259   bft_printf_flush();
260 
261   ple_coupling_mpi_intracomm_create(MPI_COMM_WORLD,
262                                     cs_glob_mpi_comm,
263                                     sat_coupling->sat_root_rank,
264                                     &(sat_coupling->comm),
265                                     local_range,
266                                     distant_range);
267 
268   bft_printf(_("[ok]\n"));
269   bft_printf(_("  Local ranks = [%d..%d], distant ranks = [%d..%d].\n\n"),
270              local_range[0], local_range[1] - 1,
271              distant_range[0], distant_range[1] - 1);
272   bft_printf_flush();
273 
274   sat_coupling->n_sat_ranks = distant_range[1] - distant_range[0];
275   sat_coupling->sat_root_rank = distant_range[0];
276 
277 #endif
278 }
279 
280 #if defined(HAVE_MPI)
281 
282 /*----------------------------------------------------------------------------
283  * Add a Code_Saturne coupling using MPI.
284  *
285  * parameters:
286  *   builder_id    <-- Code_Saturne application id in coupling builder
287  *   sat_root_rank <-- root rank associated with Code_Saturne
288  *   n_sat_ranks   <-- number of ranks associated with Code_Saturne
289  *----------------------------------------------------------------------------*/
290 
291 static void
_sat_add_mpi(int builder_id,int sat_root_rank,int n_sat_ranks)292 _sat_add_mpi(int builder_id,
293              int sat_root_rank,
294              int n_sat_ranks)
295 {
296   cs_sat_coupling_t *sat_coupling = NULL;
297   _cs_sat_coupling_builder_t *scb = _sat_coupling_builder + builder_id;
298 
299   /* Similarly to SYRTHES, we might be able to add
300      Code_Saturne couplings directly (without resorting
301      to a temporary builder), then match communications */
302 
303   cs_sat_coupling_add(scb->face_cpl_sel_c,
304                       scb->cell_cpl_sel_c,
305                       scb->face_loc_sel_c,
306                       scb->cell_loc_sel_c,
307                       scb->app_name,
308                       scb->verbosity);
309 
310   sat_coupling = cs_sat_coupling_by_id(cs_sat_coupling_n_couplings() - 1);
311 
312   sat_coupling->sat_root_rank = sat_root_rank;
313   sat_coupling->n_sat_ranks = n_sat_ranks;
314 
315   _init_comm(sat_coupling,
316              builder_id);
317 }
318 
319 /*----------------------------------------------------------------------------
320  * Print information on identified Code_Saturne couplings using MPI.
321  *
322  * This function requires coupling_builder information, and must thus
323  * be called before removing matched builder entries.
324  *----------------------------------------------------------------------------*/
325 
326 static void
_print_all_mpi_sat(void)327 _print_all_mpi_sat(void)
328 {
329   int i;
330 
331   const ple_coupling_mpi_set_t *mpi_apps = cs_coupling_get_mpi_apps();
332   const char empty_string[] = "";
333 
334   /* Loop on defined Code_Saturne instances */
335 
336   for (i = 0; i < _sat_coupling_builder_size; i++) {
337 
338     _cs_sat_coupling_builder_t *scb = _sat_coupling_builder + i;
339 
340     if (scb->match_id > -1) {
341 
342       const char *local_name = empty_string;
343       const char *distant_name = empty_string;
344 
345       const ple_coupling_mpi_set_info_t
346         ai = ple_coupling_mpi_set_get_info(mpi_apps, scb->match_id);
347 
348       if (scb->app_name != NULL)
349         local_name = scb->app_name;
350       if (ai.app_name != NULL)
351         distant_name = ai.app_name;
352 
353       bft_printf(_(" Code_Saturne coupling:\n"
354                    "   coupling id:              %d\n"
355                    "   local name:               \"%s\"\n"
356                    "   distant application name: \"%s\"\n"
357                    "   MPI application id:       %d\n"
358                    "   MPI root rank:            %d\n"
359                    "   number of MPI ranks:      %d\n\n"),
360                  i, local_name, distant_name,
361                  scb->match_id, ai.root_rank, ai.n_ranks);
362     }
363   }
364 
365   bft_printf_flush();
366 }
367 
368 /*----------------------------------------------------------------------------
369  * Initialize MPI Code_Saturne couplings using MPI.
370  *
371  * This function may be called once all couplings have been defined,
372  * and it will match defined couplings with available applications.
373  *----------------------------------------------------------------------------*/
374 
375 static void
_init_all_mpi_sat(void)376 _init_all_mpi_sat(void)
377 {
378   int i;
379 
380   int n_apps = 0;
381   int n_matched_apps = 0;
382   int n_sat_apps = 0;
383 
384   const ple_coupling_mpi_set_t *mpi_apps = cs_coupling_get_mpi_apps();
385 
386   if (mpi_apps == NULL)
387     return;
388 
389   n_apps = ple_coupling_mpi_set_n_apps(mpi_apps);
390 
391   /* First pass to count available Code_Saturne couplings */
392 
393   for (i = 0; i < n_apps; i++) {
394     const ple_coupling_mpi_set_info_t
395       ai = ple_coupling_mpi_set_get_info(mpi_apps, i);
396     if (strncmp(ai.app_type, "Code_Saturne", 12) == 0)
397       n_sat_apps += 1;
398   }
399 
400   /* In single-coupling mode, no identification necessary */
401 
402   if (n_sat_apps == 2 && _sat_coupling_builder_size == 1) {
403 
404     const int local_app_id = ple_coupling_mpi_set_get_app_id(mpi_apps);
405 
406     for (i = 0; i < n_apps; i++) {
407       const ple_coupling_mpi_set_info_t
408         ai = ple_coupling_mpi_set_get_info(mpi_apps, i);
409       if (   strncmp(ai.app_type, "Code_Saturne", 12) == 0
410           && i != local_app_id) {
411         _sat_coupling_builder->match_id = i;
412       }
413     }
414 
415     n_matched_apps += 1;
416 
417   }
418 
419   /* In multiple-coupling mode, identification is necessary */
420 
421   else {
422 
423     int j;
424     ple_coupling_mpi_set_info_t ai;
425 
426     int *sat_appinfo = NULL;
427 
428     /* First, build an array of matched/unmatched Code_Saturne applications,
429        with 2 entries per instance: matched indicator, app_id */
430 
431     BFT_MALLOC(sat_appinfo, n_sat_apps*2, int);
432 
433     n_sat_apps = 0;
434 
435     for (i = 0; i < n_apps; i++) {
436       ai = ple_coupling_mpi_set_get_info(mpi_apps, i);
437       if (strncmp(ai.app_type, "Code_Saturne", 12) == 0) {
438         sat_appinfo[n_sat_apps*2] = 0;
439         sat_appinfo[n_sat_apps*2 + 1] = i;
440         n_sat_apps += 1;
441       }
442     }
443 
444     /* Loop on defined Code_Saturne instances */
445 
446     for (i = 0; i < _sat_coupling_builder_size; i++) {
447 
448       _cs_sat_coupling_builder_t *scb = _sat_coupling_builder + i;
449 
450       /* Loop on available Code_Saturne instances to match app_names */
451 
452       if (scb->app_name != NULL) {
453 
454         for (j = 0; j < n_sat_apps; j++) {
455 
456           if (sat_appinfo[j*2] != 0) /* Consider only unmatched applications */
457             continue;
458 
459           ai = ple_coupling_mpi_set_get_info(mpi_apps, sat_appinfo[j*2 + 1]);
460           if (ai.app_name == NULL)
461             continue;
462 
463           if (strcmp(ai.app_name, scb->app_name) == 0) {
464             scb->match_id = sat_appinfo[j*2 + 1];
465             sat_appinfo[j*2] = i;
466             n_matched_apps += 1;
467             break;
468           }
469 
470         }
471 
472       }
473 
474     } /* End of loop on defined Code_Saturne instances */
475 
476     BFT_FREE(sat_appinfo);
477 
478   } /* End of test on single or multiple Code_Saturne matching algorithm */
479 
480   /* Print matching info */
481 
482   _print_all_mpi_sat();
483 
484   /* Now initialize matched couplings */
485   /*----------------------------------*/
486 
487   for (i = 0; i < _sat_coupling_builder_size; i++) {
488 
489     _cs_sat_coupling_builder_t *scb = _sat_coupling_builder + i;
490 
491     if (scb->match_id > -1) {
492       const ple_coupling_mpi_set_info_t
493         ai = ple_coupling_mpi_set_get_info(mpi_apps, scb->match_id);
494 
495       if (strncmp(ai.app_type, "Code_Saturne", 12) == 0)
496         _sat_add_mpi(i, ai.root_rank, ai.n_ranks);
497     }
498 
499   }
500 
501   /* Cleanup */
502 
503   _remove_matched_builder_entries();
504 }
505 
506 #endif /* defined(HAVE_MPI) */
507 
508 /*----------------------------------------------------------------------------
509  * Destroy a coupling structure
510  *
511  * parameters:
512  *   couplage <-> pointer to coupling structure to destroy
513  *
514  * returns:
515  *   NULL pointer
516  *----------------------------------------------------------------------------*/
517 
518 static cs_sat_coupling_t *
_sat_coupling_destroy(cs_sat_coupling_t * couplage)519 _sat_coupling_destroy(cs_sat_coupling_t  *couplage)
520 {
521   BFT_FREE(couplage->sat_name);
522 
523   BFT_FREE(couplage->face_cpl_sel);
524   BFT_FREE(couplage->cell_cpl_sel);
525   BFT_FREE(couplage->face_loc_sel);
526   BFT_FREE(couplage->cell_loc_sel);
527 
528   ple_locator_destroy(couplage->localis_cel);
529   ple_locator_destroy(couplage->localis_fbr);
530 
531   if (couplage->cells_sup != NULL)
532     fvm_nodal_destroy(couplage->cells_sup);
533   if (couplage->faces_sup != NULL)
534     fvm_nodal_destroy(couplage->faces_sup);
535 
536   BFT_FREE(couplage->distant_dist_fbr);
537   BFT_FREE(couplage->distant_of);
538   BFT_FREE(couplage->local_of);
539   BFT_FREE(couplage->distant_pond_fbr);
540   BFT_FREE(couplage->local_pond_fbr);
541 
542 #if defined(HAVE_MPI)
543   if (   couplage->comm != MPI_COMM_WORLD
544       && couplage->comm != cs_glob_mpi_comm)
545     MPI_Comm_free(&(couplage->comm));
546 #endif
547 
548   BFT_FREE(couplage);
549 
550   return NULL;
551 }
552 
553 /*----------------------------------------------------------------------------
554  * Computed some quantities needed for a centered-like interpolation
555  *  - distance JJ' for distant boundary faces
556  *  - local weighting coefficients
557  *----------------------------------------------------------------------------*/
558 
559 static void
_sat_coupling_interpolate(cs_sat_coupling_t * couplage)560 _sat_coupling_interpolate(cs_sat_coupling_t  *couplage)
561 {
562   int    icoo;
563   int    reverse;
564 
565   cs_lnum_t    ind;
566   cs_lnum_t    iel;
567   cs_lnum_t    ifac;
568 
569   cs_lnum_t    n_fbr_loc  = 0;
570   cs_lnum_t    n_fbr_dist = 0;
571 
572   cs_real_t   pdt_scal;
573   cs_real_t   surface;
574 
575   cs_real_t   distance_fbr_cel;
576   cs_real_t   distance_cel_cel;
577 
578   cs_real_t   dist_cel_fbr[3];
579   cs_real_t   vect_surf_norm[3];
580 
581   cs_real_t  *local_surf     = NULL;
582   cs_real_t  *local_xyzcen   = NULL;
583   cs_real_t  *distant_surf   = NULL;
584   cs_real_t  *distant_xyzcen = NULL;
585 
586   const cs_lnum_t   *lstfbr        = NULL;
587   const cs_lnum_t   *element       = NULL;
588   const cs_coord_t  *distant_coord = NULL;
589 
590   const cs_mesh_t  *mesh = cs_glob_mesh;
591   const cs_mesh_quantities_t  *mesh_quantities = cs_glob_mesh_quantities;
592 
593   /* Removing the connectivity and localization informations in case of
594      coupling update */
595 
596   if (couplage->distant_dist_fbr != NULL)
597     BFT_FREE(couplage->distant_dist_fbr);
598   if (couplage->distant_of != NULL)
599     BFT_FREE(couplage->distant_of);
600   if (couplage->local_of != NULL)
601     BFT_FREE(couplage->local_of);
602   if (couplage->distant_pond_fbr != NULL)
603     BFT_FREE(couplage->distant_pond_fbr);
604   if (couplage->local_pond_fbr != NULL)
605     BFT_FREE(couplage->local_pond_fbr);
606 
607   /* Interpolation structure */
608 
609   n_fbr_loc  = ple_locator_get_n_interior(couplage->localis_fbr);
610   lstfbr     = ple_locator_get_interior_list(couplage->localis_fbr);
611 
612   n_fbr_dist    = ple_locator_get_n_dist_points(couplage->localis_fbr);
613   element       = ple_locator_get_dist_locations(couplage->localis_fbr);
614   distant_coord = ple_locator_get_dist_coords(couplage->localis_fbr);
615 
616 
617   /* Calculation of the distance DJJPB defining the distance from */
618   /* the local cell center to the distant boundary face norm      */
619   /*--------------------------------------------------------------*/
620 
621   BFT_MALLOC(couplage->distant_dist_fbr, 3*n_fbr_dist, cs_real_t);
622 
623   /* Store the local surface vector of the coupled boundary faces */
624 
625   BFT_MALLOC(local_surf, 3*n_fbr_loc, cs_real_t);
626 
627   for (ind = 0 ; ind < n_fbr_loc ; ind++) {
628 
629     ifac = lstfbr[ind] - 1;
630 
631     for (icoo = 0 ; icoo < 3 ; icoo++)
632       local_surf[ind*3 + icoo] = mesh_quantities->b_face_normal[ifac*3 + icoo];
633 
634  }
635 
636   /* Get the distant faces surface vector (reverse = 1) */
637 
638   reverse = 1;
639 
640   BFT_MALLOC(distant_surf, 3*n_fbr_dist, cs_real_t);
641 
642   ple_locator_exchange_point_var(couplage->localis_fbr,
643                                  distant_surf,
644                                  local_surf,
645                                  NULL,
646                                  sizeof(cs_real_t),
647                                  3,
648                                  reverse);
649 
650   BFT_FREE(local_surf);
651 
652   /* Calculation of the JJ' vectors */
653 
654   BFT_MALLOC(distant_xyzcen, 3*n_fbr_dist, cs_real_t);
655 
656   for (ind = 0; ind < n_fbr_dist; ind++) {
657 
658     iel = element[ind] - 1;
659 
660     surface = 0.;
661     for (icoo = 0; icoo < 3; icoo++)
662       surface += distant_surf[ind*3 + icoo]*distant_surf[ind*3 + icoo];
663     surface = sqrt(surface);
664 
665     pdt_scal = 0.;
666     for (icoo = 0; icoo < 3; icoo++) {
667 
668       dist_cel_fbr[icoo] =
669         distant_coord[ind*3 + icoo] - mesh_quantities->cell_cen[iel*3 + icoo];
670 
671       /* Store the distant coordinates to compute the weighting coefficients */
672       distant_xyzcen[ind*3 + icoo] = mesh_quantities->cell_cen[iel*3 + icoo];
673 
674       vect_surf_norm[icoo] =
675         distant_surf[ind*3 + icoo] / surface;
676 
677       pdt_scal += dist_cel_fbr[icoo]*vect_surf_norm[icoo];
678 
679     }
680 
681     for (icoo = 0; icoo < 3; icoo++)
682       couplage->distant_dist_fbr[ind*3 + icoo] =
683         dist_cel_fbr[icoo] - pdt_scal*vect_surf_norm[icoo];
684 
685   }
686 
687   BFT_FREE(distant_surf);
688 
689 
690   /* Calculation of the local weighting coefficient */
691   /*------------------------------------------------*/
692 
693   BFT_MALLOC(couplage->distant_pond_fbr, n_fbr_dist, cs_real_t);
694   BFT_MALLOC(couplage->local_pond_fbr, n_fbr_loc, cs_real_t);
695 
696   /* Get the cell center coordinates (reverse = 0) */
697 
698   reverse = 0;
699 
700   BFT_MALLOC(local_xyzcen, 3*n_fbr_loc, cs_real_t);
701 
702   ple_locator_exchange_point_var(couplage->localis_fbr,
703                                  distant_xyzcen,
704                                  local_xyzcen,
705                                  NULL,
706                                  sizeof(cs_real_t),
707                                  3,
708                                  reverse);
709 
710   BFT_FREE(distant_xyzcen);
711 
712   /* Calculation of the local weighting coefficients */
713 
714   for (ind = 0 ; ind < n_fbr_loc ; ind++) {
715 
716     ifac = lstfbr[ind] - 1;
717     iel  = mesh->b_face_cells[ifac];
718 
719     surface = 0.;
720 
721     distance_fbr_cel = 0.;
722     distance_cel_cel = 0.;
723 
724     for (icoo = 0 ; icoo < 3 ; icoo++) {
725 
726       surface += mesh_quantities->b_face_normal[ifac*3 + icoo]
727         * mesh_quantities->b_face_normal[ifac*3 + icoo];
728 
729       distance_fbr_cel += mesh_quantities->b_face_normal[ifac*3 + icoo] *
730         (local_xyzcen[ind*3 + icoo] - mesh_quantities->b_face_cog[ifac*3 + icoo]);
731 
732       distance_cel_cel += mesh_quantities->b_face_normal[ifac*3 + icoo] *
733         (local_xyzcen[ind*3 + icoo] - mesh_quantities->cell_cen[iel*3 + icoo]);
734 
735     }
736 
737     surface = sqrt(surface);
738 
739     distance_fbr_cel /= surface;
740     distance_cel_cel /= surface;
741 
742     if (fabs(distance_cel_cel) > 1.e-12)
743       couplage->local_pond_fbr[ind] = distance_fbr_cel / distance_cel_cel;
744     else
745       couplage->local_pond_fbr[ind] = 0.5;
746 
747   }
748 
749   /* Get the distant weighting coefficients (reverse = 1) */
750 
751   reverse = 1;
752 
753   ple_locator_exchange_point_var(couplage->localis_fbr,
754                                  couplage->distant_pond_fbr,
755                                  couplage->local_pond_fbr,
756                                  NULL,
757                                  sizeof(cs_real_t),
758                                  1,
759                                  reverse);
760 
761   /* Calculation of the OF distance */
762   /*--------------------------------*/
763 
764   BFT_MALLOC(couplage->distant_of, 3*n_fbr_dist, cs_real_t);
765   BFT_MALLOC(couplage->local_of, 3*n_fbr_loc, cs_real_t);
766 
767   for (ind = 0 ; ind < n_fbr_loc ; ind++) {
768 
769     ifac = lstfbr[ind] - 1;
770     iel  = mesh->b_face_cells[ifac];
771 
772     surface = 0.;
773 
774     distance_fbr_cel = 0.;
775     distance_cel_cel = 0.;
776 
777     for (icoo = 0 ; icoo < 3 ; icoo++) {
778 
779       surface += mesh_quantities->b_face_normal[ifac*3 + icoo]
780         * mesh_quantities->b_face_normal[ifac*3 + icoo];
781 
782       distance_fbr_cel += mesh_quantities->b_face_normal[ifac*3 + icoo] *
783         (local_xyzcen[ind*3+icoo] - mesh_quantities->b_face_cog[ifac*3+icoo]);
784 
785       distance_cel_cel += mesh_quantities->b_face_normal[ifac*3 + icoo] *
786         (local_xyzcen[ind*3+icoo] - mesh_quantities->cell_cen[iel*3+icoo]);
787 
788     }
789 
790     surface = sqrt(surface);
791 
792     distance_fbr_cel /= surface;
793     distance_cel_cel /= surface;
794 
795     const cs_real_3_t *restrict b_face_cog
796       = (const cs_real_3_t *restrict)mesh_quantities->b_face_cog;
797     const cs_real_3_t *restrict b_face_normal
798       = (const cs_real_3_t *restrict)mesh_quantities->b_face_normal;
799 
800     for (icoo = 0 ; icoo < 3 ; icoo++) {
801       couplage->local_of[ind*3 + icoo]
802         =     b_face_cog[ifac][icoo]
803           -  (      b_face_cog[ifac][icoo] /*  O'  */
804               +     b_face_normal[ifac][icoo] *distance_fbr_cel/surface    /*J'=F+n*FJ'*/
805               - 0.5*b_face_normal[ifac][icoo] *distance_cel_cel/surface);  /*-n*I'J'/2*/
806     }
807   }
808 
809   reverse = 1;
810 
811   ple_locator_exchange_point_var(couplage->localis_fbr,
812                                  couplage->distant_of,
813                                  couplage->local_of,
814                                  NULL,
815                                  sizeof(cs_real_t),
816                                  3,
817                                  reverse);
818 
819   BFT_FREE(local_xyzcen);
820 }
821 
822 /*! (DOXYGEN_SHOULD_SKIP_THIS) \endcond */
823 
824 /*============================================================================
825  * Public function definitions for Fortran API
826  *============================================================================*/
827 
828 /*----------------------------------------------------------------------------
829  * Get number of code couplings
830  *
831  * Fortran interface:
832  *
833  * SUBROUTINE NBCCPL
834  * *****************
835  *
836  * INTEGER          NBRCPL         : <-- : number of code couplings
837  *----------------------------------------------------------------------------*/
838 
CS_PROCF(nbccpl,NBCCPL)839 void CS_PROCF (nbccpl, NBCCPL)
840 (
841  int   *n_couplings
842 )
843 {
844   if (_cs_glob_n_sat_cp < 0) {
845     _cs_glob_n_sat_cp = cs_sat_coupling_n_couplings();
846     if (_sat_coupling_builder_size > 0)
847       _cs_glob_n_sat_cp += _sat_coupling_builder_size;
848   }
849 
850   *n_couplings = _cs_glob_n_sat_cp;
851 }
852 
853 /*----------------------------------------------------------------------------
854  * Set the list of cells and boundary faces associated to a coupling
855  * and a cloud of point.
856  *
857  * The local "support" cells and boundary faces are used to localize
858  * the values in the distant "coupled" cells and faces.
859  * Depending on the role of sender and/or receiver of the current process
860  * in the coupling, some of these sets can be empty or not.
861  *
862  * The cell values are always localized and interpolated on the distant
863  * "cells" support. The face values are localized and interpolated on
864  * the distant "face" support if present, or on the distant "cell" support
865  * if not.
866  *
867  * If the input arrays LCESUP and LFBSUP are not ordered, they will be
868  * orderd in output.
869  *
870  * Fortran interface:
871  *
872  * SUBROUTINE DEFLOC
873  * *****************
874  *
875  * INTEGER          NUMCPL         : --> : coupling number
876  *----------------------------------------------------------------------------*/
877 
CS_PROCF(defloc,DEFLOC)878 void CS_PROCF (defloc, DEFLOC)
879 (
880  const int  *numcpl
881 )
882 {
883   cs_lnum_t  ind;
884   cs_lnum_t  nbr_fbr_cpl = 0, nbr_cel_cpl = 0;
885 
886   int  indic_glob[2] = {0, 0};
887   int  indic_loc[2] = {0, 0};
888 
889   char coupled_mesh_name[64];
890   cs_lnum_t *c_elt_list = NULL;
891   cs_lnum_t *f_elt_list = NULL;
892   cs_sat_coupling_t  *coupl = NULL;
893   fvm_nodal_t  *support_fbr = NULL;
894   cs_mesh_quantities_t  *mesh_quantities = cs_glob_mesh_quantities;
895 
896   int locator_options[PLE_LOCATOR_N_OPTIONS];
897   locator_options[PLE_LOCATOR_NUMBERING] = 1;
898 
899   int *point_tag = NULL;
900 
901   /* Initializations and verifications */
902 
903   if (*numcpl < 1 || *numcpl > cs_glob_sat_n_couplings)
904     bft_error(__FILE__, __LINE__, 0,
905               _("Impossible coupling number %d; there are %d couplings"),
906               *numcpl, cs_glob_sat_n_couplings);
907   else
908     coupl = cs_glob_sat_couplings[*numcpl - 1];
909 
910   /* Removing the connectivity and localization informations in case of
911      coupling update */
912 
913   if (coupl->cells_sup != NULL) fvm_nodal_destroy(coupl->cells_sup);
914   if (coupl->faces_sup != NULL) fvm_nodal_destroy(coupl->faces_sup);
915 
916   /* Create the local lists */
917 
918   if (coupl->cell_loc_sel != NULL) {
919 
920     BFT_MALLOC(c_elt_list, cs_glob_mesh->n_cells, cs_lnum_t);
921 
922     cs_selector_get_cell_num_list(coupl->cell_loc_sel,
923                                   &(coupl->nbr_cel_sup),
924                                   c_elt_list);
925 
926   }
927 
928   if (coupl->face_loc_sel != NULL) {
929 
930     BFT_MALLOC(f_elt_list, cs_glob_mesh->n_b_faces, cs_lnum_t);
931 
932     cs_selector_get_b_face_num_list(coupl->face_loc_sel,
933                                     &(coupl->nbr_fbr_sup),
934                                     f_elt_list);
935 
936   }
937 
938   if (coupl->nbr_cel_sup > 0) indic_loc[0] = 1;
939   if (coupl->nbr_fbr_sup > 0) indic_loc[1] = 1;
940 
941   for (ind = 0 ; ind < 2 ; ind++)
942     indic_glob[ind] = indic_loc[ind];
943 
944 #if defined(HAVE_MPI)
945   if (cs_glob_n_ranks > 1)
946     MPI_Allreduce (indic_loc, indic_glob, 2, MPI_INT, MPI_MAX,
947                    cs_glob_mpi_comm);
948 #endif
949 
950   if (indic_glob[0] > 0) {
951 
952     sprintf(coupled_mesh_name, _("coupled_cells_%d"), *numcpl);
953 
954     coupl->cells_sup = cs_mesh_connect_cells_to_nodal(cs_glob_mesh,
955                                                       coupled_mesh_name,
956                                                       false,
957                                                       coupl->nbr_cel_sup,
958                                                       c_elt_list);
959 
960   }
961 
962   if (indic_glob[1] > 0) {
963 
964     sprintf(coupled_mesh_name, _("coupled_faces_%d"), *numcpl);
965 
966     coupl->faces_sup = cs_mesh_connect_faces_to_nodal(cs_glob_mesh,
967                                                       coupled_mesh_name,
968                                                       false,
969                                                       0,
970                                                       coupl->nbr_fbr_sup,
971                                                       NULL,
972                                                       f_elt_list);
973 
974   }
975 
976   if (coupl->cell_loc_sel != NULL) BFT_FREE(c_elt_list);
977   if (coupl->face_loc_sel != NULL) BFT_FREE(f_elt_list);
978 
979   /* Build and initialize associated locator */
980 
981 #if defined(PLE_HAVE_MPI)
982 
983   if (coupl->localis_cel == NULL)
984     coupl->localis_cel = ple_locator_create(coupl->comm,
985                                             coupl->n_sat_ranks,
986                                             coupl->sat_root_rank);
987 
988   if (coupl->localis_fbr == NULL)
989     coupl->localis_fbr = ple_locator_create(coupl->comm,
990                                             coupl->n_sat_ranks,
991                                             coupl->sat_root_rank);
992 
993 #else
994 
995   if (coupl->localis_cel == NULL)
996     coupl->localis_cel = ple_locator_create();
997 
998   if (coupl->localis_fbr == NULL)
999     coupl->localis_fbr = ple_locator_create();
1000 
1001 #endif
1002 
1003   /* Initialization of the distant point localization */
1004 
1005   if (coupl->cell_cpl_sel != NULL) {
1006 
1007     BFT_MALLOC(c_elt_list, cs_glob_mesh->n_cells, cs_lnum_t);
1008 
1009     cs_selector_get_cell_num_list(coupl->cell_cpl_sel,
1010                                   &nbr_cel_cpl,
1011                                   c_elt_list);
1012 
1013   }
1014 
1015   if (coupl->tag_func != NULL) {
1016     BFT_MALLOC(point_tag, nbr_cel_cpl, int);
1017     coupl->tag_func(coupl->tag_context,
1018                     coupl->cells_sup,
1019                     nbr_cel_cpl,
1020                     1,
1021                     c_elt_list,
1022                     point_tag);
1023   }
1024 
1025   ple_locator_set_mesh(coupl->localis_cel,
1026                        coupl->cells_sup,
1027                        locator_options,
1028                        0.,
1029                        coupl->tolerance,
1030                        3,
1031                        nbr_cel_cpl,
1032                        c_elt_list,
1033                        point_tag,
1034                        mesh_quantities->cell_cen,
1035                        NULL,
1036                        cs_coupling_mesh_extents,
1037                        cs_coupling_point_in_mesh_p);
1038 
1039   BFT_FREE(point_tag);
1040 
1041   if (coupl->cell_cpl_sel != NULL) BFT_FREE(c_elt_list);
1042 
1043   if (coupl->face_cpl_sel != NULL) {
1044 
1045     BFT_MALLOC(f_elt_list, cs_glob_mesh->n_b_faces, cs_lnum_t);
1046 
1047     cs_selector_get_b_face_num_list(coupl->face_cpl_sel,
1048                                     &nbr_fbr_cpl,
1049                                     f_elt_list);
1050 
1051   }
1052 
1053   if (indic_glob[1] > 0)
1054     support_fbr = coupl->faces_sup;
1055   else
1056     support_fbr = coupl->cells_sup;
1057 
1058   if (coupl->tag_func != NULL) {
1059     BFT_MALLOC(point_tag, nbr_fbr_cpl, int);
1060     coupl->tag_func(coupl->tag_context,
1061                     support_fbr,
1062                     nbr_fbr_cpl,
1063                     1,
1064                     f_elt_list,
1065                     point_tag);
1066   }
1067 
1068   ple_locator_set_mesh(coupl->localis_fbr,
1069                        support_fbr,
1070                        locator_options,
1071                        0.,
1072                        coupl->tolerance,
1073                        3,
1074                        nbr_fbr_cpl,
1075                        f_elt_list,
1076                        point_tag,
1077                        mesh_quantities->b_face_cog,
1078                        NULL,
1079                        cs_coupling_mesh_extents,
1080                        cs_coupling_point_in_mesh_p);
1081 
1082   BFT_FREE(point_tag);
1083 
1084   if (coupl->face_cpl_sel != NULL) BFT_FREE(f_elt_list);
1085 
1086   /* Computed some quantities needed for a centered-like interpolation */
1087 
1088   if (coupl->localis_fbr != NULL)
1089     _sat_coupling_interpolate(coupl);
1090 
1091 #if 0
1092   /* TODO: associate the FVM meshes to the post-processing,
1093      with a fonction giving a pointer to the associated FVM structures,
1094      and another enabling its compacting or removing */
1095   {
1096     fvm_writer_t *w = fvm_writer_init("coupled_mesh",
1097                                       NULL,
1098                                       "EnSight Gold",
1099                                       "binary",
1100                                       FVM_WRITER_FIXED_MESH);
1101 
1102     fvm_writer_export_nodal(w, coupl->cells_sup);
1103     fvm_writer_finalize(w);
1104 
1105   }
1106 #endif
1107 
1108   /* Compacting the interpolation support (could be removed) */
1109 
1110   if (coupl->cells_sup != NULL)
1111     fvm_nodal_reduce(coupl->cells_sup, 1);
1112   if (coupl->faces_sup != NULL)
1113     fvm_nodal_reduce(coupl->faces_sup, 1);
1114 
1115 #if 0 && defined(DEBUG) && !defined(NDEBUG)
1116   ple_locator_dump(coupl->localis_cel);
1117   ple_locator_dump(coupl->localis_fbr);
1118 #endif
1119 }
1120 
1121 /*----------------------------------------------------------------------------
1122  * Get the number of cells and boundary faces, "support", coupled and not
1123  * localized associated to a given coupling
1124  *
1125  * Fortran interface:
1126  *
1127  * SUBROUTINE NBECPL
1128  * *****************
1129  *
1130  * INTEGER          NUMCPL         : --> : coupling number
1131  * INTEGER          NCESUP         : <-- : number of "support" cells
1132  * INTEGER          NFBSUP         : <-- : number of "support" boundary faces
1133  * INTEGER          NCECPL         : <-- : number of coupled cells
1134  * INTEGER          NFBCPL         : <-- : number of coupled boundary faces
1135  * INTEGER          NCENCP         : <-- : number of not coupled cells
1136  *                                 :     : (since not localized)
1137  * INTEGER          NFBNCP         : <-- : number of not coupled boundary faces
1138  *                                 :     : (since not localized)
1139  *----------------------------------------------------------------------------*/
1140 
CS_PROCF(nbecpl,NBECPL)1141 void CS_PROCF (nbecpl, NBECPL)
1142 (
1143  const int        *numcpl,
1144        cs_lnum_t  *ncesup,
1145        cs_lnum_t  *nfbsup,
1146        cs_lnum_t  *ncecpl,
1147        cs_lnum_t  *nfbcpl,
1148        cs_lnum_t  *ncencp,
1149        cs_lnum_t  *nfbncp
1150 )
1151 {
1152   cs_sat_coupling_t  *coupl = NULL;
1153 
1154   /* Initializations and verifications */
1155 
1156   if (*numcpl < 1 || *numcpl > cs_glob_sat_n_couplings)
1157     bft_error(__FILE__, __LINE__, 0,
1158               _("Impossible coupling number %d; there are %d couplings"),
1159               *numcpl, cs_glob_sat_n_couplings);
1160   else
1161     coupl = cs_glob_sat_couplings[*numcpl - 1];
1162 
1163   *ncesup = coupl->nbr_cel_sup;
1164   *nfbsup = coupl->nbr_fbr_sup;
1165 
1166   *ncecpl = 0;
1167   *nfbcpl = 0;
1168 
1169   *ncencp = 0;
1170   *nfbncp = 0;
1171 
1172   if (coupl->localis_cel != NULL) {
1173     *ncecpl = ple_locator_get_n_interior(coupl->localis_cel);
1174     *ncencp = ple_locator_get_n_exterior(coupl->localis_cel);
1175   }
1176 
1177   if (coupl->localis_fbr != NULL) {
1178     *nfbcpl = ple_locator_get_n_interior(coupl->localis_fbr);
1179     *nfbncp = ple_locator_get_n_exterior(coupl->localis_fbr);
1180   }
1181 
1182 }
1183 
1184 /*----------------------------------------------------------------------------
1185  * Get the lists of coupled cells and boundary faces (i.e. receiving)
1186  * associated to a given coupling
1187  *
1188  * The number of cells and boundary faces, got with NBECPL(), are used
1189  * for arguments coherency checks.
1190  *
1191  * Fortran interface:
1192  *
1193  * SUBROUTINE LELCPL
1194  * *****************
1195  *
1196  * INTEGER          NUMCPL         : --> : coupling number
1197  * INTEGER          NCECPL         : --> : number of coupled cells
1198  * INTEGER          NFBCPL         : --> : number of coupled boundary faces
1199  * INTEGER          LCECPL(*)      : <-- : list of coupled cells
1200  * INTEGER          LFBCPL(*)      : <-- : list of coupled boundary faces
1201  *----------------------------------------------------------------------------*/
1202 
CS_PROCF(lelcpl,LELCPL)1203 void CS_PROCF (lelcpl, LELCPL)
1204 (
1205  const int        *numcpl,
1206  const cs_lnum_t  *ncecpl,
1207  const cs_lnum_t  *nfbcpl,
1208        cs_lnum_t  *lcecpl,
1209        cs_lnum_t  *lfbcpl
1210 )
1211 {
1212   cs_lnum_t  ind;
1213 
1214   cs_lnum_t  _ncecpl = 0;
1215   cs_lnum_t  _nfbcpl = 0;
1216 
1217   cs_sat_coupling_t  *coupl = NULL;
1218 
1219   const cs_lnum_t  *lst = NULL;
1220 
1221   /* Initializations and verifications */
1222 
1223   if (*numcpl < 1 || *numcpl > cs_glob_sat_n_couplings)
1224     bft_error(__FILE__, __LINE__, 0,
1225               _("Impossible coupling number %d; there are %d couplings"),
1226               *numcpl, cs_glob_sat_n_couplings);
1227   else
1228     coupl = cs_glob_sat_couplings[*numcpl - 1];
1229 
1230   if (coupl->localis_cel != NULL)
1231     _ncecpl = ple_locator_get_n_interior(coupl->localis_cel);
1232 
1233   if (coupl->localis_fbr != NULL)
1234     _nfbcpl = ple_locator_get_n_interior(coupl->localis_fbr);
1235 
1236   if (*ncecpl != _ncecpl || *nfbcpl != _nfbcpl)
1237     bft_error(__FILE__, __LINE__, 0,
1238               _("Coupling %d: inconsistent arguments for LELCPL()\n"
1239                 "NCECPL = %d and NFBCPL = %d are indicated.\n"
1240                 "The values for this coupling should be %d and %d."),
1241               *numcpl, (int)(*ncecpl), (int)(*nfbcpl),
1242               (int)_ncecpl, (int)_nfbcpl);
1243 
1244   /* Copy lists (would be useless with a pure C API) */
1245 
1246   if (_ncecpl > 0) {
1247     lst = ple_locator_get_interior_list(coupl->localis_cel);
1248     for (ind = 0 ; ind < _ncecpl ; ind++)
1249       lcecpl[ind] = lst[ind];
1250   }
1251 
1252   if (_nfbcpl > 0) {
1253     lst = ple_locator_get_interior_list(coupl->localis_fbr);
1254     for (ind = 0 ; ind < _nfbcpl ; ind++)
1255       lfbcpl[ind] = lst[ind];
1256   }
1257 }
1258 
1259 /*----------------------------------------------------------------------------
1260  * Get the lists of not coupled cells and boundary faces (i.e. receiving but
1261  * not localized) associated to a given coupling
1262  *
1263  * The number of cells and boundary faces, got with NBECPL(), are used
1264  * for arguments coherency checks.
1265  *
1266  * Fortran interface:
1267  *
1268  * SUBROUTINE LENCPL
1269  * *****************
1270  *
1271  * INTEGER          NUMCPL         : --> : coupling number
1272  * INTEGER          NCENCP         : --> : number of not coupled cells
1273  * INTEGER          NFBNCP         : --> : number of not coupled boundary faces
1274  * INTEGER          LCENCP(*)      : <-- : list of not coupled cells
1275  * INTEGER          LFBNCP(*)      : <-- : list of not coupled boundary faces
1276  *----------------------------------------------------------------------------*/
1277 
CS_PROCF(lencpl,LENCPL)1278 void CS_PROCF (lencpl, LENCPL)
1279 (
1280  const int        *numcpl,
1281  const cs_lnum_t  *ncencp,
1282  const cs_lnum_t  *nfbncp,
1283        cs_lnum_t  *lcencp,
1284        cs_lnum_t  *lfbncp
1285 )
1286 {
1287   cs_lnum_t  ind;
1288 
1289   cs_lnum_t  _ncencp = 0;
1290   cs_lnum_t  _nfbncp = 0;
1291   cs_sat_coupling_t  *coupl = NULL;
1292 
1293   const cs_lnum_t  *lst = NULL;
1294 
1295   /* Initializations and verifications */
1296 
1297   if (*numcpl < 1 || *numcpl > cs_glob_sat_n_couplings)
1298     bft_error(__FILE__, __LINE__, 0,
1299               _("Impossible coupling number %d; there are %d couplings"),
1300               *numcpl, cs_glob_sat_n_couplings);
1301   else
1302     coupl = cs_glob_sat_couplings[*numcpl - 1];
1303 
1304   if (coupl->localis_cel != NULL)
1305     _ncencp = ple_locator_get_n_exterior(coupl->localis_cel);
1306 
1307   if (coupl->localis_fbr != NULL)
1308     _nfbncp = ple_locator_get_n_exterior(coupl->localis_fbr);
1309 
1310   if (*ncencp != _ncencp || *nfbncp != _nfbncp)
1311     bft_error(__FILE__, __LINE__, 0,
1312               _("Coupling %d: inconsistent arguments for LELNCP()\n"
1313                 "NCENCP = %d and NFBNCP = %d are indicated.\n"
1314                 "The values for this coupling should be %d and %d."),
1315               *numcpl, (int)(*ncencp), (int)(*nfbncp),
1316               (int)_ncencp, (int)_nfbncp);
1317 
1318   /* Copy lists (would be useless with a pure C API) */
1319 
1320   if (_ncencp > 0) {
1321     lst = ple_locator_get_exterior_list(coupl->localis_cel);
1322     for (ind = 0 ; ind < _ncencp ; ind++)
1323       lcencp[ind] = lst[ind];
1324   }
1325 
1326   if (_nfbncp > 0) {
1327     lst = ple_locator_get_exterior_list(coupl->localis_fbr);
1328     for (ind = 0 ; ind < _nfbncp ; ind++)
1329       lfbncp[ind] = lst[ind];
1330   }
1331 }
1332 
1333 /*----------------------------------------------------------------------------
1334  * Get the number of distant point associated to a given coupling
1335  * and localized on the local domain
1336  *
1337  * Fortran interface:
1338  *
1339  * SUBROUTINE NPDCPL
1340  * *****************
1341  *
1342  * INTEGER          NUMCPL         : --> : coupling number
1343  * INTEGER          NCEDIS         : <-- : number of distant cells
1344  * INTEGER          NFBDIS         : <-- : numbre de distant boundary faces
1345  *----------------------------------------------------------------------------*/
1346 
CS_PROCF(npdcpl,NPDCPL)1347 void CS_PROCF (npdcpl, NPDCPL)
1348 (
1349  const int        *numcpl,
1350        cs_lnum_t  *ncedis,
1351        cs_lnum_t  *nfbdis
1352 )
1353 {
1354   cs_sat_coupling_t  *coupl = NULL;
1355 
1356   /* Verifications */
1357 
1358   if (*numcpl < 1 || *numcpl > cs_glob_sat_n_couplings)
1359     bft_error(__FILE__, __LINE__, 0,
1360               _("Impossible coupling number %d; there are %d couplings"),
1361               *numcpl, cs_glob_sat_n_couplings);
1362   else
1363     coupl = cs_glob_sat_couplings[*numcpl - 1];
1364 
1365   /* Get the number of points */
1366 
1367   *ncedis = 0;
1368   *nfbdis = 0;
1369 
1370   if (coupl->localis_cel != NULL)
1371     *ncedis = ple_locator_get_n_dist_points(coupl->localis_cel);
1372 
1373   if (coupl->localis_fbr != NULL)
1374     *nfbdis = ple_locator_get_n_dist_points(coupl->localis_fbr);
1375 
1376 }
1377 
1378 /*----------------------------------------------------------------------------
1379  * Get the distant points coordinates associated to a given coupling
1380  * and a list of points, and the elements number and type (cell or face)
1381  * "containing" this points.
1382  *
1383  * The number of distant points NBRPTS must be equal to one the arguments
1384  * NCEDIS or NFBDIS given by NPDCPL(), and is given here for coherency checks
1385  * between the arguments NUMCPL and ITYSUP.
1386  *
1387  * Fortran interface:
1388  *
1389  * SUBROUTINE COOCPL
1390  * *****************
1391  *
1392  * INTEGER          NUMCPL         : --> : coupling number
1393  * INTEGER          NBRPTS         : --> : number of distant points
1394  * INTEGER          ITYDIS         : --> : 1 : access to the points associated
1395  *                                 :     :     to the distant cells
1396  *                                 :     : 2 : access to the points associated
1397  *                                 :     :     to the distant boundary faces
1398  * INTEGER          ITYLOC         : <-- : 1 : localization on the local cells
1399  *                                 :     : 2 : localization on the local faces
1400  * INTEGER          LOCPTS(*)      : <-- : "containing" number associated to
1401  *                                 :     :   each point
1402  * DOUBLE PRECISION COOPTS(3,*)    : <-- : distant point coordinates
1403  * DOUBLE PRECISION DJPPTS(3,*)    : <-- : distant vectors to the coupled face
1404  * DOUBLE PRECISION PNDPTS(*)      : <-- : distant weighting coefficients
1405  *----------------------------------------------------------------------------*/
1406 
CS_PROCF(coocpl,COOCPL)1407 void CS_PROCF (coocpl, COOCPL)
1408 (
1409  const int        *numcpl,
1410  const cs_lnum_t  *nbrpts,
1411  const int        *itydis,
1412        int        *ityloc,
1413        cs_lnum_t  *locpts,
1414        cs_real_t  *coopts,
1415        cs_real_t  *djppts,
1416        cs_real_t  *dofpts,
1417        cs_real_t  *pndpts
1418 )
1419 {
1420   cs_lnum_t  ind, icoo;
1421 
1422   cs_lnum_t  n_pts_dist = 0;
1423   cs_sat_coupling_t  *coupl = NULL;
1424   ple_locator_t  *localis = NULL;
1425 
1426   /* Initializations and verifications */
1427 
1428   if (*numcpl < 1 || *numcpl > cs_glob_sat_n_couplings)
1429     bft_error(__FILE__, __LINE__, 0,
1430               _("Impossible coupling number %d; there are %d couplings"),
1431               *numcpl, cs_glob_sat_n_couplings);
1432   else
1433     coupl = cs_glob_sat_couplings[*numcpl - 1];
1434 
1435   *ityloc = 0;
1436 
1437   if (*itydis == 1) {
1438     localis = coupl->localis_cel;
1439     *ityloc = 1;
1440   }
1441   else if (*itydis == 2) {
1442     localis = coupl->localis_fbr;
1443     if (coupl->nbr_fbr_sup > 0)
1444       *ityloc = 2;
1445     else
1446       *ityloc = 1;
1447   }
1448 
1449   if (localis != NULL)
1450     n_pts_dist = ple_locator_get_n_dist_points(localis);
1451 
1452   if (*nbrpts != n_pts_dist)
1453     bft_error(__FILE__, __LINE__, 0,
1454               _("Coupling %d: inconsistent arguments for COOCPL()\n"
1455                 "ITYDIS = %d and NBRPTS = %d are indicated.\n"
1456                 "The value for NBRPTS should be %d."),
1457               *numcpl, (int)(*itydis), (int)(*nbrpts), (int)n_pts_dist);
1458 
1459   /* Creation the local lists */
1460 
1461   if (localis != NULL) {
1462 
1463     n_pts_dist = ple_locator_get_n_dist_points(localis);
1464 
1465     if (n_pts_dist > 0) {
1466 
1467       const cs_lnum_t   *element;
1468       const cs_coord_t  *coord;
1469 
1470       element = ple_locator_get_dist_locations(localis);
1471       coord   = ple_locator_get_dist_coords(localis);
1472 
1473       for (ind = 0 ; ind < n_pts_dist ; ind++) {
1474         locpts[ind] = element[ind];
1475         for (icoo = 0 ; icoo < 3 ; icoo++)
1476           coopts[ind*3 + icoo] = coord[ind*3 + icoo];
1477       }
1478 
1479       if (*itydis == 2)
1480         for (ind = 0 ; ind < n_pts_dist ; ind++)
1481           for (icoo = 0 ; icoo < 3 ; icoo++) {
1482             djppts[ind*3 + icoo] = coupl->distant_dist_fbr[ind*3 + icoo];
1483             dofpts[ind*3 + icoo] = coupl->distant_of[ind*3 + icoo];
1484             pndpts[ind] = coupl->distant_pond_fbr[ind];
1485           }
1486 
1487     }
1488 
1489   }
1490 
1491 }
1492 
1493 /*----------------------------------------------------------------------------
1494  * Get the weighting coefficient needed for a centered-like interpolation
1495  * in the case of a coupling on boundary faces.
1496  *
1497  * Fortran interface:
1498  *
1499  * SUBROUTINE PONDCP
1500  * *****************
1501  *
1502  * INTEGER          NUMCPL         : --> : coupling number
1503  * INTEGER          NBRPTS         : --> : number of distant points
1504  * INTEGER          ITYLOC         : <-- : 1 : localization on the local cells
1505  *                                 :     : 2 : localization on the local faces
1506  * DOUBLE PRECISION PNDCPL(*)      : <-- : weighting coefficients
1507  *----------------------------------------------------------------------------*/
1508 
CS_PROCF(pondcp,PONDCP)1509 void CS_PROCF (pondcp, PONDCP)
1510 (
1511  const int        *numcpl,
1512  const cs_lnum_t  *nbrpts,
1513        int        *ityloc,
1514        cs_real_t  *pndcpl,
1515        cs_real_t  *distof
1516 )
1517 {
1518   cs_lnum_t       nfbcpl = 0;
1519   cs_sat_coupling_t  *coupl = NULL;
1520   ple_locator_t  *localis = NULL;
1521 
1522   /* Initializations and verifications */
1523 
1524   if (*numcpl < 1 || *numcpl > cs_glob_sat_n_couplings)
1525     bft_error(__FILE__, __LINE__, 0,
1526               _("Impossible coupling number %d; there are %d couplings"),
1527               *numcpl, cs_glob_sat_n_couplings);
1528   else
1529     coupl = cs_glob_sat_couplings[*numcpl - 1];
1530 
1531   if (*ityloc == 1)
1532     bft_error(__FILE__, __LINE__, 0,
1533               _("The centered interpolation scheme is not available\n"
1534                 "when coupling cells"));
1535   else if (*ityloc == 2)
1536     localis = coupl->localis_fbr;
1537 
1538 
1539   if (localis != NULL)
1540     nfbcpl = ple_locator_get_n_interior(localis);
1541 
1542   if (*nbrpts != nfbcpl)
1543     bft_error(__FILE__, __LINE__, 0,
1544               _("Coupling %d: inconsistent arguments for PNDCPL().\n"
1545                 "ITYLOC = %d and NBRPTS = %d are indicated.\n"
1546                 "NBRPTS should be %d."),
1547               *numcpl, (int)(*ityloc), (int)(*nbrpts), (int)nfbcpl);
1548 
1549   /* Creation of the local lists */
1550 
1551   if (localis != NULL) {
1552 
1553     if (nfbcpl > 0) {
1554 
1555       for (cs_lnum_t ind = 0 ; ind < nfbcpl ; ind++) {
1556         pndcpl[ind] = coupl->local_pond_fbr[ind];
1557         for (cs_lnum_t icoo = 0 ; icoo < 3 ; icoo++)
1558           distof[ind*3 + icoo] = coupl->local_of[ind*3 + icoo];
1559       }
1560 
1561     }
1562 
1563   }
1564 }
1565 
1566 /*----------------------------------------------------------------------------
1567  * Exchange a variable associated to a set of point and a coupling.
1568  *
1569  * Fortran interface:
1570  *
1571  * SUBROUTINE VARCPL
1572  * *****************
1573  *
1574  * INTEGER          NUMCPL         : --> : coupling number
1575  * INTEGER          NBRDIS         : --> : number of values to send
1576  * INTEGER          NBRLOC         : --> : number of values to receive
1577  * INTEGER          ITYVAR         : --> : 1 : variables defined at cells
1578  *                                 :     : 2 : variables defined at faces
1579  * INTEGER          STRIDE         : --> : 1 : for scalars
1580  *                                 :     : 3 : for vectors
1581  * DOUBLE PRECISION VARDIS(*)      : --> : distant variable(to send)
1582  * DOUBLE PRECISION VARLOC(*)      : <-- : local variable (to receive)
1583  *----------------------------------------------------------------------------*/
1584 
CS_PROCF(varcpl,VARCPL)1585 void CS_PROCF (varcpl, VARCPL)
1586 (
1587  const int        *numcpl,
1588  const cs_lnum_t  *nbrdis,
1589  const cs_lnum_t  *nbrloc,
1590  const int        *ityvar,
1591  const cs_lnum_t  *stride,
1592        cs_real_t  *vardis,
1593        cs_real_t  *varloc
1594 )
1595 {
1596   cs_lnum_t  n_val_dist_ref = 0;
1597   cs_lnum_t  n_val_loc_ref = 0;
1598   cs_real_t  *val_dist = NULL;
1599   cs_real_t  *val_loc = NULL;
1600   cs_sat_coupling_t  *coupl = NULL;
1601   ple_locator_t  *localis = NULL;
1602 
1603   /* Initializations and verifications */
1604 
1605   if (*numcpl < 1 || *numcpl > cs_glob_sat_n_couplings)
1606     bft_error(__FILE__, __LINE__, 0,
1607               _("Impossible coupling number %d; there are %d couplings"),
1608               *numcpl, cs_glob_sat_n_couplings);
1609   else
1610     coupl = cs_glob_sat_couplings[*numcpl - 1];
1611 
1612   if (*ityvar == 1)
1613     localis = coupl->localis_cel;
1614   else if (*ityvar == 2)
1615     localis = coupl->localis_fbr;
1616 
1617   if (localis != NULL) {
1618     n_val_dist_ref = ple_locator_get_n_dist_points(localis);
1619     n_val_loc_ref  = ple_locator_get_n_interior(localis);
1620   }
1621 
1622   if (*nbrdis > 0 && *nbrdis != n_val_dist_ref)
1623     bft_error(__FILE__, __LINE__, 0,
1624               _("Coupling %d: inconsistent arguments for VARCPL()\n"
1625                 "ITYVAR = %d and NBRDIS = %d are indicated.\n"
1626                 "NBRDIS should be 0 or %d."),
1627               *numcpl, (int)(*ityvar), (int)(*nbrdis), (int)n_val_dist_ref);
1628 
1629   if (*nbrloc > 0 && *nbrloc != n_val_loc_ref)
1630     bft_error(__FILE__, __LINE__, 0,
1631               _("Coupling %d: inconsistent arguments for VARCPL()\n"
1632                 "ITYVAR = %d and NBRLOC = %d are indicated.\n"
1633                 "NBRLOC should be 0 or %d."),
1634               *numcpl, (int)(*ityvar), (int)(*nbrloc), (int)n_val_loc_ref);
1635 
1636   /* Create the local lists */
1637 
1638   if (localis != NULL) {
1639 
1640     if (*nbrdis > 0)
1641       val_dist = vardis;
1642     if (*nbrloc > 0)
1643       val_loc = varloc;
1644 
1645     ple_locator_exchange_point_var(localis,
1646                                    val_dist,
1647                                    val_loc,
1648                                    NULL,
1649                                    sizeof(cs_real_t),
1650                                    *stride,
1651                                    0);
1652   }
1653 
1654 }
1655 
1656 /*----------------------------------------------------------------------------
1657  * Array of integers exchange, associated to a given coupling.
1658  *
1659  * It is assumed that the arrays have the same size and the same values on
1660  * each group of processus (local and distant).
1661  *
1662  * Fortran interface:
1663  *
1664  * SUBROUTINE TBICPL
1665  * *****************
1666  *
1667  * INTEGER          NUMCPL         : --> : coupling number
1668  * INTEGER          NBRDIS         : --> : number of values to send
1669  * INTEGER          NBRLOC         : --> : number of values to receive
1670  * INTEGER          TABDIS(*)      : --> : distant values (to send)
1671  * INTEGER          TABLOC(*)      : <-- : local values (to receive)
1672  *----------------------------------------------------------------------------*/
1673 
CS_PROCF(tbicpl,TBICPL)1674 void CS_PROCF (tbicpl, TBICPL)
1675 (
1676  const int        *numcpl,
1677  const cs_lnum_t  *nbrdis,
1678  const cs_lnum_t  *nbrloc,
1679        cs_lnum_t  *vardis,
1680        cs_lnum_t  *varloc
1681 )
1682 {
1683   cs_lnum_t  nbr = 0;
1684   bool  distant = false;
1685 
1686 #if defined(HAVE_MPI)
1687 
1688   MPI_Status  status;
1689 
1690   cs_sat_coupling_t  *coupl = NULL;
1691 
1692   /* Initializations and verifications */
1693 
1694   if (*numcpl < 1 || *numcpl > cs_glob_sat_n_couplings)
1695     bft_error(__FILE__, __LINE__, 0,
1696               _("Impossible coupling number %d; there are %d couplings"),
1697               *numcpl, cs_glob_sat_n_couplings);
1698   else
1699     coupl = cs_glob_sat_couplings[*numcpl - 1];
1700 
1701   if (coupl->comm != MPI_COMM_NULL) {
1702 
1703     distant = true;
1704 
1705     /* Exchange between the groups master node */
1706 
1707     if (cs_glob_rank_id < 1)
1708       MPI_Sendrecv(vardis, *nbrdis, CS_MPI_LNUM, coupl->sat_root_rank, 0,
1709                    varloc, *nbrloc, CS_MPI_LNUM, coupl->sat_root_rank, 0,
1710                    coupl->comm, &status);
1711 
1712     /* Synchronization inside a group */
1713 
1714     if (cs_glob_n_ranks > 1)
1715       MPI_Bcast (varloc, *nbrloc, CS_MPI_LNUM, 0, cs_glob_mpi_comm);
1716 
1717   }
1718 
1719 #endif /* defined(HAVE_MPI) */
1720 
1721   if (distant == false) {
1722 
1723     nbr = CS_MIN(*nbrdis, *nbrloc);
1724 
1725     for (cs_lnum_t ind = 0; ind < nbr; ind++)
1726       varloc[ind] = vardis[ind];
1727 
1728   }
1729 }
1730 
1731 /*----------------------------------------------------------------------------
1732  * Array of reals exchange, associated to a given coupling.
1733  *
1734  * It is assumed that the arrays have the same size and the same values on
1735  * each group of processus (local and distant).
1736  *
1737  * Fortran interface:
1738  *
1739  * SUBROUTINE TBRCPL
1740  * *****************
1741  *
1742  * INTEGER          NUMCPL         : --> : coupling number
1743  * INTEGER          NBRDIS         : --> : number of values to send
1744  * INTEGER          NBRLOC         : --> : number of values to receive
1745  * DOUBLE PRECISION TABDIS(*)      : --> : distant values (to send)
1746  * DOUBLE PRECISION TABLOC(*)      : <-- : local values (to receive)
1747  *----------------------------------------------------------------------------*/
1748 
CS_PROCF(tbrcpl,TBRCPL)1749 void CS_PROCF (tbrcpl, TBRCPL)
1750 (
1751  const int        *numcpl,
1752  const cs_lnum_t  *nbrdis,
1753  const cs_lnum_t  *nbrloc,
1754        cs_real_t  *vardis,
1755        cs_real_t  *varloc
1756 )
1757 {
1758   cs_lnum_t  nbr = 0;
1759   bool  distant = false;
1760 
1761 #if defined(HAVE_MPI)
1762 
1763   MPI_Status  status;
1764 
1765   cs_sat_coupling_t  *coupl = NULL;
1766 
1767   /* Initializations and verifications */
1768 
1769   if (*numcpl < 1 || *numcpl > cs_glob_sat_n_couplings)
1770     bft_error(__FILE__, __LINE__, 0,
1771               _("Impossible coupling number %d; there are %d couplings"),
1772               *numcpl, cs_glob_sat_n_couplings);
1773   else
1774     coupl = cs_glob_sat_couplings[*numcpl - 1];
1775 
1776   if (coupl->comm != MPI_COMM_NULL) {
1777 
1778     distant = true;
1779 
1780     /* Exchange between the groups master node */
1781 
1782     if (cs_glob_rank_id < 1)
1783       MPI_Sendrecv(vardis, *nbrdis, CS_MPI_REAL, coupl->sat_root_rank, 0,
1784                    varloc, *nbrloc, CS_MPI_REAL, coupl->sat_root_rank, 0,
1785                    coupl->comm, &status);
1786 
1787     /* Synchronization inside a group */
1788 
1789     if (cs_glob_n_ranks > 1)
1790       MPI_Bcast(varloc, *nbrloc, CS_MPI_REAL, 0, cs_glob_mpi_comm);
1791 
1792   }
1793 
1794 #endif /* defined(HAVE_MPI) */
1795 
1796   if (distant == false) {
1797 
1798     nbr = CS_MIN(*nbrdis, *nbrloc);
1799 
1800     for (cs_lnum_t ind = 0; ind < nbr; ind++)
1801       varloc[ind] = vardis[ind];
1802 
1803   }
1804 }
1805 
1806 /*----------------------------------------------------------------------------
1807  * Compute the maximum value of an integer variable associated to a coupling.
1808  *
1809  * It is assumed that the integer value is the same for each group of
1810  * processus (local and distant).
1811  *
1812  * Fortran interface:
1813  *
1814  * SUBROUTINE MXICPL
1815  * *****************
1816  *
1817  * INTEGER          NUMCPL         : --> : coupling number
1818  * INTEGER          VALDIS         : --> : distant value (to send)
1819  * INTEGER          VALMAX         : <-- : local maximum (to receive)
1820  *----------------------------------------------------------------------------*/
1821 
CS_PROCF(mxicpl,MXICPL)1822 void CS_PROCF (mxicpl, MXICPL)
1823 (
1824  const int        *const numcpl,
1825        cs_lnum_t  *const vardis,
1826        cs_lnum_t  *const varmax
1827 )
1828 {
1829   bool  distant = false;
1830 
1831 #if defined(HAVE_MPI)
1832 
1833   cs_sat_coupling_t  *coupl = NULL;
1834 
1835   /* Initializations and verifications */
1836 
1837   if (*numcpl < 1 || *numcpl > cs_glob_sat_n_couplings)
1838     bft_error(__FILE__, __LINE__, 0,
1839               _("Impossible coupling number %d; there are %d couplings"),
1840               *numcpl, cs_glob_sat_n_couplings);
1841   else
1842     coupl = cs_glob_sat_couplings[*numcpl - 1];
1843 
1844   if (coupl->comm != MPI_COMM_NULL) {
1845 
1846     distant = true;
1847 
1848     MPI_Allreduce(vardis, varmax, 1, CS_MPI_LNUM, MPI_MAX, coupl->comm);
1849 
1850   }
1851 
1852 #endif /* defined(HAVE_MPI) */
1853 
1854   if (distant == false) {
1855 
1856     *varmax = *vardis;
1857 
1858   }
1859 }
1860 
1861 /*============================================================================
1862  * Public function definitions
1863  *============================================================================*/
1864 
1865 /*----------------------------------------------------------------------------*/
1866 /*!
1867  * \brief Define new Code_Saturne coupling.
1868  *
1869  * The arguments to \ref cs_sat_coupling_define are:
1870  * \param[in] saturne_name          matching Code_Saturne application name
1871  * \param[in] boundary_cpl_criteria boundary face selection criteria for coupled
1872  *                                  faces, or NULL
1873  * \param[in] volume_cpl_criteria   cell selection criteria for coupled cells, or
1874                                     NULL
1875  * \param[in] boundary_loc_criteria boundary face selection criteria for location
1876  *                                  (not functional)
1877  * \param[in] volume_loc_criteria   cell selection criteria for location
1878  * \param[in] verbosity             verbosity level
1879  *
1880  * In the case of only 2 Code_Saturne instances, the 'saturne_name' argument
1881  * is ignored, as there is only one matching possibility.
1882  *
1883  * In case of multiple couplings, a coupling will be matched with available
1884  * Code_Saturne instances based on the 'saturne_name' argument.
1885  */
1886 /*----------------------------------------------------------------------------*/
1887 
1888 void
cs_sat_coupling_define(const char * saturne_name,const char * boundary_cpl_criteria,const char * volume_cpl_criteria,const char * boundary_loc_criteria,const char * volume_loc_criteria,int verbosity)1889 cs_sat_coupling_define(const char  *saturne_name,
1890                        const char  *boundary_cpl_criteria,
1891                        const char  *volume_cpl_criteria,
1892                        const char  *boundary_loc_criteria,
1893                        const char  *volume_loc_criteria,
1894                        int          verbosity)
1895 {
1896   _cs_sat_coupling_builder_t *scb = NULL;
1897 
1898   /* Add corresponding coupling to temporary Code_Saturne couplings array */
1899 
1900   BFT_REALLOC(_sat_coupling_builder,
1901               _sat_coupling_builder_size + 1,
1902               _cs_sat_coupling_builder_t);
1903 
1904   scb = &(_sat_coupling_builder[_sat_coupling_builder_size]);
1905 
1906   scb->match_id = -1;
1907 
1908   scb->app_name = NULL;
1909   if (saturne_name != NULL) {
1910     BFT_MALLOC(scb->app_name, strlen(saturne_name) + 1, char);
1911     strcpy(scb->app_name, saturne_name);
1912   }
1913 
1914   scb->face_cpl_sel_c = NULL;
1915   if (boundary_cpl_criteria != NULL) {
1916     BFT_MALLOC(scb->face_cpl_sel_c, strlen(boundary_cpl_criteria) + 1, char);
1917     strcpy(scb->face_cpl_sel_c, boundary_cpl_criteria);
1918   }
1919 
1920   scb->cell_cpl_sel_c = NULL;
1921   if (volume_cpl_criteria != NULL) {
1922     BFT_MALLOC(scb->cell_cpl_sel_c, strlen(volume_cpl_criteria) + 1, char);
1923     strcpy(scb->cell_cpl_sel_c, volume_cpl_criteria);
1924   }
1925 
1926   scb->face_loc_sel_c = NULL;
1927   if (boundary_loc_criteria != NULL) {
1928     BFT_MALLOC(scb->face_loc_sel_c, strlen(boundary_loc_criteria) + 1, char);
1929     strcpy(scb->face_loc_sel_c, boundary_loc_criteria);
1930   }
1931 
1932   scb->cell_loc_sel_c = NULL;
1933   if (volume_loc_criteria != NULL) {
1934     BFT_MALLOC(scb->cell_loc_sel_c, strlen(volume_loc_criteria) + 1, char);
1935     strcpy(scb->cell_loc_sel_c, volume_loc_criteria);
1936   }
1937 
1938   scb->verbosity = verbosity;
1939 
1940   _sat_coupling_builder_size += 1;
1941 }
1942 
1943 /*----------------------------------------------------------------------------
1944  * Get number of Code_Saturne couplings.
1945  *
1946  * returns:
1947  *   number of Code_Saturne couplings
1948  *----------------------------------------------------------------------------*/
1949 
1950 int
cs_sat_coupling_n_couplings(void)1951 cs_sat_coupling_n_couplings(void)
1952 {
1953   return cs_glob_sat_n_couplings;
1954 }
1955 
1956 /*----------------------------------------------------------------------------
1957  * Get pointer to Code_Saturne coupling.
1958  *
1959  * parameters:
1960  *   coupling_id <-- Id (0 to n-1) of Code_Saturne coupling
1961  *
1962  * returns:
1963  *   pointer to Code_Saturne coupling structure
1964  *----------------------------------------------------------------------------*/
1965 
1966 cs_sat_coupling_t *
cs_sat_coupling_by_id(int coupling_id)1967 cs_sat_coupling_by_id(int coupling_id)
1968 {
1969   cs_sat_coupling_t  *retval = NULL;
1970 
1971   if (   coupling_id > -1
1972       && coupling_id < cs_glob_sat_n_couplings)
1973     retval = cs_glob_sat_couplings[coupling_id];
1974 
1975   return retval;
1976 }
1977 
1978 /*----------------------------------------------------------------------------
1979  * Initialize Code_Saturne couplings.
1980  *
1981  * This function may be called once all couplings have been defined,
1982  * and it will match defined couplings with available applications.
1983  *----------------------------------------------------------------------------*/
1984 
1985 void
cs_sat_coupling_all_init(void)1986 cs_sat_coupling_all_init(void)
1987 {
1988   /* First, try using other MPI communicators */
1989 
1990 #if defined(HAVE_MPI)
1991 
1992   if (_sat_coupling_builder_size > 0)
1993     _init_all_mpi_sat();
1994 
1995 #endif
1996 
1997   /* Print unmatched instances */
1998 
1999   if (_sat_coupling_builder_size > 0) {
2000 
2001     bft_printf("Unmatched Code_Saturne couplings:\n"
2002                "---------------------------------\n\n");
2003 
2004     _print_all_unmatched_sat();
2005 
2006     bft_error(__FILE__, __LINE__, 0,
2007               _("At least 1 Code_Saturne coupling was defined for which\n"
2008                 "no communication with a Code_Saturne instance is possible."));
2009   }
2010 }
2011 
2012 /*----------------------------------------------------------------------------
2013  * Create a sat_coupling_t structure.
2014  *
2015  * parameters:
2016  *   ref_axis           <-- reference axis
2017  *   face_sel_criterion <-- criterion for selection of boundary faces
2018  *   cell_sel_criterion <-- criterion for selection of cells
2019  *   sat_name           <-- Code_Saturne application name
2020  *   verbosity          <-- verbosity level
2021  *----------------------------------------------------------------------------*/
2022 
2023 void
cs_sat_coupling_add(const char * face_cpl_sel_c,const char * cell_cpl_sel_c,const char * face_loc_sel_c,const char * cell_loc_sel_c,const char * sat_name,int verbosity)2024 cs_sat_coupling_add(const char  *face_cpl_sel_c,
2025                     const char  *cell_cpl_sel_c,
2026                     const char  *face_loc_sel_c,
2027                     const char  *cell_loc_sel_c,
2028                     const char  *sat_name,
2029                     int          verbosity)
2030 {
2031   cs_sat_coupling_t *sat_coupling = NULL;
2032 
2033   /* Allocate _cs_sat_coupling_t structure */
2034 
2035   BFT_REALLOC(cs_glob_sat_couplings,
2036               cs_glob_sat_n_couplings + 1, cs_sat_coupling_t*);
2037   BFT_MALLOC(sat_coupling, 1, cs_sat_coupling_t);
2038 
2039   sat_coupling->sat_name = NULL;
2040   sat_coupling->tag_func = NULL;
2041   sat_coupling->tag_context = NULL;
2042 
2043   if (sat_name != NULL) {
2044     BFT_MALLOC(sat_coupling->sat_name, strlen(sat_name) + 1, char);
2045     strcpy(sat_coupling->sat_name, sat_name);
2046   }
2047 
2048   /* Selection criteria  */
2049 
2050   sat_coupling->face_cpl_sel = NULL;
2051   sat_coupling->cell_cpl_sel = NULL;
2052   sat_coupling->face_loc_sel = NULL;
2053   sat_coupling->cell_loc_sel = NULL;
2054 
2055   if (face_cpl_sel_c != NULL) {
2056     BFT_MALLOC(sat_coupling->face_cpl_sel, strlen(face_cpl_sel_c) + 1, char);
2057     strcpy(sat_coupling->face_cpl_sel, face_cpl_sel_c);
2058   }
2059   if (cell_cpl_sel_c != NULL) {
2060     BFT_MALLOC(sat_coupling->cell_cpl_sel, strlen(cell_cpl_sel_c) + 1, char);
2061     strcpy(sat_coupling->cell_cpl_sel, cell_cpl_sel_c);
2062   }
2063 
2064   if (face_loc_sel_c != NULL) {
2065     BFT_MALLOC(sat_coupling->face_loc_sel, strlen(face_loc_sel_c) + 1, char);
2066     strcpy(sat_coupling->face_loc_sel, face_loc_sel_c);
2067   }
2068   if (cell_loc_sel_c != NULL) {
2069     BFT_MALLOC(sat_coupling->cell_loc_sel, strlen(cell_loc_sel_c) + 1, char);
2070     strcpy(sat_coupling->cell_loc_sel, cell_loc_sel_c);
2071   }
2072 
2073   sat_coupling->faces_sup = NULL;
2074   sat_coupling->cells_sup = NULL;
2075 
2076   sat_coupling->localis_fbr = NULL;
2077   sat_coupling->localis_cel = NULL;
2078 
2079   sat_coupling->nbr_fbr_sup = 0;
2080   sat_coupling->nbr_cel_sup = 0;
2081 
2082   sat_coupling->tolerance = 0.1;
2083   sat_coupling->verbosity = verbosity;
2084 
2085   /* Geometric quantities arrays for interpolation */
2086 
2087   sat_coupling->distant_dist_fbr = NULL;
2088   sat_coupling->distant_of = NULL;
2089   sat_coupling->local_of = NULL;
2090   sat_coupling->distant_pond_fbr = NULL;
2091   sat_coupling->local_pond_fbr = NULL;
2092 
2093   /* Initialize communicators */
2094 
2095 #if defined(HAVE_MPI)
2096 
2097   sat_coupling->comm = MPI_COMM_NULL;
2098   sat_coupling->n_sat_ranks = 0;
2099   sat_coupling->sat_root_rank = -1;
2100 
2101 #endif
2102 
2103   /* Update coupling array and return */
2104 
2105   cs_glob_sat_couplings[cs_glob_sat_n_couplings] = sat_coupling;
2106   cs_glob_sat_n_couplings++;
2107 }
2108 
2109 /*----------------------------------------------------------------------------
2110  * Create a new internal Code_Saturne coupling.
2111  *
2112  * arguments:
2113  *   tag_func          <-- pointer to tagging function
2114  *   tag_context       <-- pointer to tagging function context
2115  *   boundary_criteria <-- boundary face selection criteria, or NULL
2116  *   volume_criteria   <-- volume cell selection criteria, or NULL
2117  *   loc_tolerance     <-- location tolerance factor (0.1 recommended)
2118  *   verbosity         <-- verbosity level
2119  *----------------------------------------------------------------------------*/
2120 
2121 void
cs_sat_coupling_add_internal(cs_sat_coupling_tag_t * tag_func,void * tag_context,const char * boundary_cpl_criteria,const char * volume_cpl_criteria,const char * boundary_loc_criteria,const char * volume_loc_criteria,float loc_tolerance,int verbosity)2122 cs_sat_coupling_add_internal(cs_sat_coupling_tag_t  *tag_func,
2123                              void                   *tag_context,
2124                              const char             *boundary_cpl_criteria,
2125                              const char             *volume_cpl_criteria,
2126                              const char             *boundary_loc_criteria,
2127                              const char             *volume_loc_criteria,
2128                              float                   loc_tolerance,
2129                              int                     verbosity)
2130 {
2131   CS_UNUSED(loc_tolerance);
2132 
2133   cs_sat_coupling_add(boundary_cpl_criteria,
2134                       volume_cpl_criteria,
2135                       boundary_loc_criteria,
2136                       volume_loc_criteria,
2137                       NULL,
2138                       verbosity);
2139 
2140   cs_sat_coupling_t *sat_coupling
2141     = cs_sat_coupling_by_id(cs_sat_coupling_n_couplings() - 1);
2142 
2143   sat_coupling->tag_func = tag_func;
2144   sat_coupling->tag_context = tag_context;
2145 
2146 #if defined(HAVE_MPI)
2147 
2148   sat_coupling->comm = cs_glob_mpi_comm;
2149   sat_coupling->sat_root_rank = 0;
2150   sat_coupling->n_sat_ranks = cs_glob_n_ranks;
2151 
2152 #endif
2153 }
2154 
2155 /*----------------------------------------------------------------------------
2156  * Destroy all couplings
2157  *----------------------------------------------------------------------------*/
2158 
2159 void
cs_sat_coupling_all_finalize(void)2160 cs_sat_coupling_all_finalize(void)
2161 {
2162   int  i;
2163 
2164   for (i = 0 ; i < cs_glob_sat_n_couplings ; i++)
2165     _sat_coupling_destroy(cs_glob_sat_couplings[i]);
2166 
2167   BFT_FREE(cs_glob_sat_couplings);
2168 
2169   cs_glob_sat_n_couplings = 0;
2170 }
2171 
2172 /*----------------------------------------------------------------------------*/
2173 
2174 END_C_DECLS
2175