1
/*============================================================================
2
* Unit test for fvm_interface.c;
3
*============================================================================*/
6
This file is part of Code_Saturne, a general-purpose CFD tool.
8
Copyright (C) 1998-2011 EDF S.A.
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
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
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.
25
/*----------------------------------------------------------------------------*/
27
#if defined(HAVE_CONFIG_H)
28
#include "cs_config.h"
37
#include <bft_error.h>
39
#include <bft_printf.h>
41
#include "fvm_interface.h"
42
#include "fvm_parall.h"
43
#include "fvm_periodicity.h"
45
/*---------------------------------------------------------------------------*/
47
static fvm_periodicity_t *
48
_init_periodicity(void)
50
double t1[3] = {1., 0., 0.};
51
double t2[3] = {0., 1., 0.};
52
double t3[3] = {0., 0., 1.};
54
/* double t4[3] = {1.02, 0., 0.}; */
55
/* double a1[3] = {0., 0., 1.}; */
56
/* double i1[3] = {0., 0., 0.}; */
58
fvm_periodicity_t *p = NULL;
60
p = fvm_periodicity_create(0.1);
62
fvm_periodicity_add_translation(p, 1, t1);
63
fvm_periodicity_add_translation(p, 2, t2);
64
fvm_periodicity_add_translation(p, 3, t3);
66
/* fvm_periodicity_add_rotation(p, 2, 90., a1, i1); */
67
/* fvm_periodicity_add_translation(p, 3, t3); */
69
fvm_periodicity_combine(p, 1);
74
static fvm_interface_set_t *
75
_periodic_is(int comm_size,
79
const fvm_periodicity_t *perio)
81
fvm_gnum_t ii, jj, kk, couple_id;
83
fvm_gnum_t n_max = n_side * n_side * n_side;
84
fvm_lnum_t n_elements = ((n_max+1)*4) / (3*comm_size);
85
int periodicity_num[3] = {1, 2, 3};
87
fvm_lnum_t *n_periodic_couples = NULL;
88
fvm_gnum_t **periodic_couples = NULL;
89
fvm_gnum_t *global_number = NULL;
91
fvm_interface_set_t *ifset = NULL;
95
BFT_MALLOC(global_number, n_elements, fvm_gnum_t);
97
for (ii = 0; ii < (fvm_gnum_t)n_elements; ii++) {
98
global_number[ii] = (n_elements * 3 / 4) * comm_rank + ii + 1;
100
if (comm_rank == comm_size -1 && global_number[n_elements - 1] < n_max)
101
global_number[n_elements - 1] = n_max;
105
BFT_MALLOC(n_periodic_couples, n_periodic_lists, fvm_lnum_t);
106
BFT_MALLOC(periodic_couples, n_periodic_lists, fvm_gnum_t *);
108
for (ii = 0; ii < (fvm_gnum_t)n_periodic_lists; ii++) {
110
n_periodic_couples[ii] = 0;
111
periodic_couples[ii] = NULL;
113
if (comm_rank == 0 || comm_rank == 1) {
115
n_periodic_couples[ii] = n_side*n_side;
116
BFT_MALLOC(periodic_couples[ii], n_periodic_couples[ii]*2, fvm_gnum_t);
122
for (jj = 0; jj < n_side*n_side; jj++) {
123
periodic_couples[ii][couple_id*2] = (jj*n_side) + 1;
124
periodic_couples[ii][couple_id*2+1]
125
= periodic_couples[ii][couple_id*2] + (n_side-1);
130
for (jj = 0; jj < n_side; jj++) {
131
for (kk = 0; kk < n_side; kk++) {
132
periodic_couples[ii][couple_id*2] = jj*n_side*n_side + kk + 1;
133
periodic_couples[ii][couple_id*2+1]
134
= periodic_couples[ii][couple_id*2] + (n_side-1)*n_side;
140
for (jj = 0; jj < n_side*n_side; jj++) {
141
periodic_couples[ii][couple_id*2] = jj + 1;
142
periodic_couples[ii][couple_id*2+1]
143
= periodic_couples[ii][couple_id*2] + (n_side-1)*n_side*n_side;
151
n_periodic_couples[ii] = couple_id;
152
BFT_REALLOC(periodic_couples[ii], n_periodic_couples[ii]*2, fvm_gnum_t);
158
ifset = fvm_interface_set_create(n_elements,
165
(const fvm_gnum_t **const)periodic_couples);
167
for (ii = 0; ii < (fvm_gnum_t)n_periodic_lists; ii++)
168
BFT_FREE(periodic_couples[ii]);
169
BFT_FREE(periodic_couples);
170
BFT_FREE(n_periodic_couples);
173
BFT_FREE(global_number);
178
/*----------------------------------------------------------------------------
179
* Fonction d'impression d'un message sur la sortie standard
180
*----------------------------------------------------------------------------*/
182
static int _bft_printf_proxy
184
const char *const format,
188
static FILE *f = NULL;
193
#if defined(HAVE_MPI)
194
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
196
sprintf (filename, "fvm_interface_test_out.%d", rank);
197
f = fopen(filename, "w");
201
return vfprintf(f, format, arg_ptr);
204
/*---------------------------------------------------------------------------*/
207
main (int argc, char *argv[])
209
char mem_trace_name[32];
212
fvm_interface_set_t *ifset = NULL;
213
fvm_periodicity_t *perio = NULL;
215
#if defined(HAVE_MPI)
223
fvm_lnum_t n_elements = 20;
224
fvm_gnum_t *global_number = NULL;
228
MPI_Init(&argc, &argv);
230
fvm_parall_set_mpi_comm(MPI_COMM_WORLD);
232
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
233
MPI_Comm_size(MPI_COMM_WORLD, &size);
236
sprintf(mem_trace_name, "fvm_interface_test_mem.%d", rank);
238
strcpy(mem_trace_name, "fvm_interface_test_mem");
239
bft_mem_init(mem_trace_name);
241
bft_printf_proxy_set(_bft_printf_proxy);
243
/* Build arbitrary interface */
245
BFT_MALLOC(global_number, n_elements, fvm_gnum_t);
247
for (ii = 0; ii < n_elements; ii++) {
248
global_number[ii] = (n_elements * 3 / 4) * rank + ii + 1;
251
ifset = fvm_interface_set_create(n_elements,
260
BFT_FREE(global_number);
262
/* Serialize dump of interfaces */
265
MPI_Recv(&sync, 1, MPI_INT, rank - 1, 0, MPI_COMM_WORLD, &status);
267
bft_printf("Interface on rank %d:\n\n", rank);
269
fvm_interface_set_dump(ifset);
272
MPI_Send(&sync, 1, MPI_INT, rank + 1, 0, MPI_COMM_WORLD);
274
/* We are finished for this interface */
276
ifset = fvm_interface_set_destroy(ifset);
280
strcpy(mem_trace_name, "fvm_interface_test_mem");
281
bft_mem_init(mem_trace_name);
283
bft_printf_proxy_set(_bft_printf_proxy);
285
#endif /* (HAVE_MPI) */
287
/* Now build interface with periodicity */
289
perio = _init_periodicity();
291
/* Dump periodicity info for first rank (identical on all ranks) */
294
fvm_periodicity_dump(perio);
296
ifset = _periodic_is(size,
298
4, /* n vertices/cube side */
299
3, /* n_periodic_lists */
302
#if defined(HAVE_MPI)
304
/* Serialize dump of interfaces */
307
MPI_Recv(&sync, 1, MPI_INT, rank - 1, 0, MPI_COMM_WORLD, &status);
309
#endif /* (HAVE_MPI) */
311
bft_printf("Periodic Interface on rank %d:\n\n", rank);
313
fvm_interface_set_dump(ifset);
315
#if defined(HAVE_MPI)
318
MPI_Send(&sync, 1, MPI_INT, rank + 1, 0, MPI_COMM_WORLD);
320
#endif /* (HAVE_MPI) */
322
/* We are finished */
324
ifset = fvm_interface_set_destroy(ifset);
325
fvm_periodicity_destroy(perio);
329
#if defined(HAVE_MPI)