1
/*============================================================================
2
* Management of the GUI parameters file: radiative transfer
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"
31
/*----------------------------------------------------------------------------
32
* Standard C library headers
33
*----------------------------------------------------------------------------*/
43
/*----------------------------------------------------------------------------
45
*----------------------------------------------------------------------------*/
48
#include <bft_error.h>
49
#include <bft_printf.h>
51
/*----------------------------------------------------------------------------
53
*----------------------------------------------------------------------------*/
55
#include "fvm_selector.h"
57
/*----------------------------------------------------------------------------
59
*----------------------------------------------------------------------------*/
62
#include "cs_gui_variables.h"
63
#include "cs_gui_util.h"
64
#include "cs_gui_boundary_conditions.h"
65
#include "cs_gui_specific_physics.h"
69
/*----------------------------------------------------------------------------
70
* Header for the current file
71
*----------------------------------------------------------------------------*/
73
#include "cs_gui_radiative_transfer.h"
75
/*----------------------------------------------------------------------------*/
79
/*=============================================================================
80
* Local Macro Definitions
81
*============================================================================*/
83
/* debugging switch */
86
/*============================================================================
87
* Local Structure Definitions
88
*============================================================================*/
90
/*----------------------------------------------------------------------------
91
* Structure associated to boundary conditions definition
92
*----------------------------------------------------------------------------*/
95
char **label; /* label for each boundary zone */
96
char **nature; /* nature for each boundary zone */
100
double *conductivity;
102
double *thermal_conductivity;
103
double *external_temp;
104
double *internal_temp;
105
double *conduction_flux;
106
} cs_radiative_boundary_t;
108
/*----------------------------------------------------------------------------
109
* Private global variables for boundary conditions
110
*----------------------------------------------------------------------------*/
112
static cs_radiative_boundary_t *boundary = NULL;
114
/*----------------------------------------------------------------------------
115
* Private global variables for the treatment
116
* of NOMVAR. NOMVAR is a characters fortran array
117
*----------------------------------------------------------------------------*/
119
static int _cs_gui_max_vars = 0;
120
static int _cs_gui_last_var = 0;
121
static char ** _cs_gui_var_rayt = NULL;
124
/*============================================================================
125
* Private function definitions
126
*============================================================================*/
128
/*-----------------------------------------------------------------------------
129
* Return integer parameters for radiation
132
* param --> name of parameter
133
* keyword <-- value of parameter
134
*----------------------------------------------------------------------------*/
137
_radiative_transfer_int(const char *const param,
143
path = cs_xpath_init_path();
144
cs_xpath_add_elements(&path, 3,
145
"thermophysical_models",
146
"radiative_transfer",
148
cs_xpath_add_function_text(&path);
150
if (cs_gui_get_int(path, &value)) *keyword = value;
155
/*-----------------------------------------------------------------------------
156
* Return float parameters for radiation
159
* param --> name of parameter
160
* keyword <-- value of parameter
161
*----------------------------------------------------------------------------*/
164
_radiative_transfer_double(const char *const param,
165
double *const keyword)
170
path = cs_xpath_init_path();
171
cs_xpath_add_elements(&path, 3,
172
"thermophysical_models",
173
"radiative_transfer",
175
cs_xpath_add_function_text(&path);
177
if (cs_gui_get_double(path, &value)) *keyword = value;
182
/*-----------------------------------------------------------------------------
183
* Return value of the parameter of the character type for radiation
186
* param --> name of parameter
187
* keyword <-- value of parameter
188
*----------------------------------------------------------------------------*/
191
_radiative_transfer_char(const char *const param,
197
path = cs_xpath_init_path();
198
cs_xpath_add_elements(&path, 3,
199
"thermophysical_models",
200
"radiative_transfer",
202
cs_xpath_add_attribute(&path, "status");
204
if(cs_gui_get_status(path, &result)) *keyword = result;
210
/*-----------------------------------------------------------------------------
211
* Return status and label of the property for post treatment of radiation
214
* name --> name of property
215
* value <-- value of status
216
*----------------------------------------------------------------------------*/
219
_radiative_transfer_char_post(const char *const name,
220
int *const list_value,
221
int *const record_value)
229
path = cs_xpath_init_path();
231
cs_xpath_add_elements(&path, 3,
232
"thermophysical_models",
233
"radiative_transfer",
235
cs_xpath_add_test_attribute(&path, "name", name);
237
BFT_MALLOC(path1, strlen(path)+1, char);
239
BFT_MALLOC(path2, strlen(path)+1, char);
242
cs_xpath_add_attribute(&path, "label");
243
label = cs_gui_get_attribute_value(path);
245
cs_xpath_add_element(&path1, "listing_printing");
246
cs_xpath_add_attribute(&path1, "status");
247
if (cs_gui_get_status(path1, &result)) {
251
cs_xpath_add_element(&path2, "postprocessing_recording");
252
cs_xpath_add_attribute(&path2, "status");
253
if (cs_gui_get_status(path2, &result)) {
264
/*-----------------------------------------------------------------------------
265
* Return value of the type of absorption coefficient for radiation
268
* param --> name of parameter "absorption coefficient"
269
* keyword <-- value of the type of the coefficent
270
*----------------------------------------------------------------------------*/
273
_radiative_transfer_type(const char *const param,
279
path = cs_xpath_init_path();
281
cs_xpath_add_elements(&path, 3,
282
"thermophysical_models",
283
"radiative_transfer",
286
cs_xpath_add_attribute(&path, "type");
288
type = cs_gui_get_attribute_value(path);
291
if (cs_gui_strcmp(type, "constant"))
293
else if (cs_gui_strcmp(type, "variable"))
295
else if (cs_gui_strcmp(type, "formula"))
297
else if (cs_gui_strcmp(type, "modak"))
300
bft_error (__FILE__, __LINE__, 0,
301
_("unknow type %s\n"), type);
308
/*----------------------------------------------------------------------------
309
* Return value of radiative variable
312
* label --> label of boundary nature
313
* param --> name of the variable
314
* value <-- value of the variable
315
*----------------------------------------------------------------------------*/
318
_radiative_boundary(const char *const label,
319
const char *const param,
325
path = cs_xpath_init_path();
326
cs_xpath_add_elements(&path, 2,
327
"boundary_conditions",
329
cs_xpath_add_test_attribute(&path, "label", label);
330
cs_xpath_add_elements(&path, 2,
333
cs_xpath_add_function_text(&path);
335
if (cs_gui_get_double(path, &res)) {
343
/*----------------------------------------------------------------------------
344
* Return int value of the type of radiative condition
347
* label --> label of boundary "wall"
348
* itpimp <-- if wall faces with imposed temperature
349
* ipgrno <-- if grey or black wall faces
350
* iprefl <-- if reflecting wall faces
351
* ifgrno <-- if grey or black wall faces and conduction flux imposed
352
* ifrefl <-- if refecting wall faces and conduction flux imposed
353
*----------------------------------------------------------------------------*/
356
_radiative_boundary_type(const char *const label,
367
path = cs_xpath_init_path();
368
cs_xpath_add_elements(&path, 2,
369
"boundary_conditions",
371
cs_xpath_add_test_attribute(&path, "label", label);
373
cs_xpath_add_element(&path, "radiative_data");
374
cs_xpath_add_attribute(&path,"choice");
375
type = cs_gui_get_attribute_value(path);
377
if (cs_gui_strcmp(type, "itpimp"))
379
else if (cs_gui_strcmp(type, "ipgrno"))
381
else if (cs_gui_strcmp(type, "iprefl"))
383
else if (cs_gui_strcmp(type, "ifgrno"))
385
else if (cs_gui_strcmp(type, "ifrefl"))
389
bft_error (__FILE__, __LINE__, 0,
390
_("Xpath request failed %s \n"), path);
398
/*----------------------------------------------------------------------------
399
* Return maximum value of output zone
400
*----------------------------------------------------------------------------*/
403
_radiative_boundary_output_zone_max(void)
406
int nb_zone, zone_max = 0;
408
path = cs_xpath_init_path();
409
cs_xpath_add_elements(&path, 4,
410
"boundary_conditions",
415
nb_zone = cs_gui_get_nb_element(path);
418
cs_xpath_add_function_text(&path);
419
zone_max = cs_gui_get_max_value(path);
427
/*-----------------------------------------------------------------------------
428
* Copy a variable name to private variable names array
431
* varname --> name or label of the variable/scalar/property
432
* ipp --> index from the fortran array associated to varname
433
*----------------------------------------------------------------------------*/
436
_cs_gui_copy_varname(const char *varname, int ipp)
440
if (ipp < 1 || ipp > _cs_gui_last_var)
441
bft_error(__FILE__, __LINE__, 0,
442
_("Variable index %d out of bounds (1 to %d)"),
443
ipp, _cs_gui_last_var);
447
if (_cs_gui_var_rayt[ipp-1] == NULL)
448
BFT_MALLOC(_cs_gui_var_rayt[ipp-1], l + 1, char);
450
else if (strlen(_cs_gui_var_rayt[ipp-1]) != l)
451
BFT_REALLOC(_cs_gui_var_rayt[ipp-1], l + 1, char);
453
strcpy(_cs_gui_var_rayt[ipp-1], varname);
456
/*============================================================================
457
* Public Fortran function definitions
458
*============================================================================*/
460
/*----------------------------------------------------------------------------
462
*----------------------------------------------------------------------------*/
464
void CS_PROCF (uiray1, UIRAY1) (int *const iirayo,
474
model = cs_gui_get_thermophysical_model("radiative_transfer");
476
if (cs_gui_strcmp(model, "off"))
478
else if (cs_gui_strcmp(model, "dom"))
480
else if (cs_gui_strcmp(model, "p-1"))
485
_radiative_transfer_char("restart", isuird);
486
_radiative_transfer_int("directions_number", ndirec);
487
_radiative_transfer_int("frequency", nfreqr);
488
_radiative_transfer_int("thermal_radiative_source_term", idiver);
489
_radiative_transfer_int("temperature_listing_printing", iimpar);
490
_radiative_transfer_int("intensity_resolution_listing_printing", iimlum);
493
bft_printf("==>UIRAY1\n");
494
bft_printf("--rayonnement : %s (iirayo = %i)\n", model, *iirayo);
497
bft_printf("--isuird = %d\n", *isuird);
498
bft_printf("--ndirec = %d\n", *ndirec);
499
bft_printf("--nfreqr = %d\n", *nfreqr);
500
bft_printf("--idiver = %d\n", *idiver);
501
bft_printf("--iimpar = %d\n", *iimpar);
502
bft_printf("--iimlum = %d\n", *iimlum);
508
/*----------------------------------------------------------------------------
510
*----------------------------------------------------------------------------*/
512
void CS_PROCF (uiray4, UIRAY4) (int *const nbrayf,
517
int list_ind, record_ind = 0;
520
const char *const _cs_properties_name2[8] = {
524
"thermal_conductivity",
532
for (i=0 ; i < *nbrayf ; i++)
536
label = _radiative_transfer_char_post(_cs_properties_name2[i], &list_ind, &record_ind);
537
irayvf[i] = record_ind;
539
_cs_gui_copy_varname(label, i + 1);
544
bft_printf("==>UIRAY4\n");
547
for (i=0 ; i < *nbrayf ; i++)
549
bft_printf(_("--output boundary faces: %s value %i \n"),
557
/*-----------------------------------------------------------------------------
558
* Indirection between the solver numbering and the XML one
559
* for physical properties of radiative transfer
560
*----------------------------------------------------------------------------*/
562
void CS_PROCF (uirapr, UIRAPR) (const int *const nprayc,
563
const int *const nprayb,
564
const int *const nrphas,
565
const int *const ipppro,
566
const int *const ipproc,
567
const int *const ilumin,
568
const int *const iqx,
569
const int *const iqy,
570
const int *const iqz,
571
const int *const itsre,
572
const int *const itsri,
573
const int *const iabs,
574
const int *const iemi,
575
const int *const icak)
582
cs_var_t *vars = cs_glob_var;
585
vars->nprop += *nprayc;
586
vars->nprayc = *nprayc;
588
BFT_REALLOC(vars->properties_ipp, vars->nprop, int);
589
BFT_REALLOC(vars->propce, vars->nprop, int);
590
BFT_REALLOC(vars->properties_name, vars->nprop, char*);
593
vars->properties_ipp[n] = ipppro[ ipproc[ *ilumin -1 ]-1 ];
594
vars->propce[n] = *ilumin;
595
BFT_MALLOC(vars->properties_name[n], strlen("intensity")+1, char);
596
strcpy(vars->properties_name[n++], "intensity");
599
vars->properties_ipp[n] = ipppro[ ipproc[ *iqx -1 ]-1 ];
600
vars->propce[n] = *iqx;
601
BFT_MALLOC(vars->properties_name[n], strlen("qrad_x")+1, char);
602
strcpy(vars->properties_name[n++], "qrad_x");
605
vars->properties_ipp[n] = ipppro[ ipproc[ *iqy -1 ]-1 ];
606
vars->propce[n] = *iqy;
607
BFT_MALLOC(vars->properties_name[n], strlen("qrad_y")+1, char);
608
strcpy(vars->properties_name[n++], "qrad_y");
611
vars->properties_ipp[n] = ipppro[ ipproc[ *iqz -1 ]-1 ];
612
vars->propce[n] = *iqz;
613
BFT_MALLOC(vars->properties_name[n], strlen("qrad_z")+1, char);
614
strcpy(vars->properties_name[n++], "qrad_z");
617
vars->properties_ipp[n] = ipppro[ ipproc[ itsre[0] -1 ]-1 ];
618
vars->propce[n] = itsre[0];
619
BFT_MALLOC(vars->properties_name[n], strlen("radiative_source_term")+1, char);
620
strcpy(vars->properties_name[n++], "radiative_source_term");
622
/* ITSRE loop on classes */
623
BFT_MALLOC(name, strlen("radiative_source_term_")+1 + 2, char);
624
BFT_MALLOC(snumpp, 1 + 2, char);
625
strcpy(name, "radiative_source_term_");
626
for (i = 1; i < *nrphas; i++)
628
sprintf(snumpp, "%2.2i", i);
629
strcat(name, snumpp);
631
vars->properties_ipp[n] = ipppro[ ipproc[ itsre[i] -1 ]-1 ];
632
vars->propce[n] = itsre[i];
633
BFT_MALLOC(vars->properties_name[n], strlen(name)+1, char);
634
strcpy(vars->properties_name[n++], name);
636
strcpy(name, "radiative_source_term_");
640
vars->properties_ipp[n] = ipppro[ ipproc[ itsri[0] -1 ]-1 ];
641
vars->propce[n] = itsri[0];
642
BFT_MALLOC(vars->properties_name[n], strlen("implicit_source_term")+1, char);
643
strcpy(vars->properties_name[n++], "implicit_source_term");
645
/* ITSRI loop on classes */
646
BFT_REALLOC(name, strlen("implicit_source_term_")+1 + 2, char);
647
strcpy(name, "implicit_source_term_");
648
for (i = 1; i < *nrphas; i++)
650
sprintf(snumpp, "%2.2i", i);
651
strcat(name, snumpp);
653
vars->properties_ipp[n] = ipppro[ ipproc[ itsri[i] -1 ]-1 ];
654
vars->propce[n] = itsri[i];
655
BFT_MALLOC(vars->properties_name[n], strlen(name)+1, char);
656
strcpy(vars->properties_name[n++], name);
658
strcpy(name, "implicit_source_term_");
662
vars->properties_ipp[n] = ipppro[ ipproc[ iabs[0] -1 ]-1 ];
663
vars->propce[n] = iabs[0];
664
BFT_MALLOC(vars->properties_name[n], strlen("absorption")+1, char);
665
strcpy(vars->properties_name[n++], "absorption");
667
/* IABS loop on classes */
668
BFT_REALLOC(name, strlen("absorption_")+1 + 2, char);
669
strcpy(name, "absorption_");
670
for (i = 1; i < *nrphas; i++)
672
sprintf(snumpp, "%2.2i", i);
673
strcat(name, snumpp);
675
vars->properties_ipp[n] = ipppro[ ipproc[ iabs[i] -1 ]-1 ];
676
vars->propce[n] = iabs[i];
677
BFT_MALLOC(vars->properties_name[n], strlen(name)+1, char);
678
strcpy(vars->properties_name[n++], name);
680
strcpy(name, "absorption_");
684
vars->properties_ipp[n] = ipppro[ ipproc[ iemi[0] -1 ]-1 ];
685
vars->propce[n] = iemi[0];
686
BFT_MALLOC(vars->properties_name[n], strlen("emission")+1, char);
687
strcpy(vars->properties_name[n++], "emission");
689
/* IEMI loop on classes */
690
BFT_REALLOC(name, strlen("emission_")+1 + 2, char);
691
strcpy(name, "emission_");
692
for (i = 1; i < *nrphas; i++)
694
sprintf(snumpp, "%2.2i", i);
695
strcat(name, snumpp);
697
vars->properties_ipp[n] = ipppro[ ipproc[ iemi[i] -1 ]-1 ];
698
vars->propce[n] = iemi[i];
699
BFT_MALLOC(vars->properties_name[n], strlen(name)+1, char);
700
strcpy(vars->properties_name[n++], name);
702
strcpy(name, "emission_");
706
vars->properties_ipp[n] = ipppro[ ipproc[ icak[0] -1 ]-1 ];
707
vars->propce[n] = icak[0];
708
BFT_MALLOC(vars->properties_name[n], strlen("absorption_coefficient")+1, char);
709
strcpy(vars->properties_name[n++], "absorption_coefficient");
711
/* ICAK loop on classes */
712
BFT_REALLOC(name, strlen("absorption_coefficient_")+1 + 2, char);
713
strcpy(name, "absorption_coefficient_");
714
for (i = 1; i < *nrphas; i++)
716
sprintf(snumpp, "%2.2i", i);
717
strcat(name, snumpp);
719
vars->properties_ipp[n] = ipppro[ ipproc[ icak[i] -1 ]-1 ];
720
vars->propce[n] = icak[i];
721
BFT_MALLOC(vars->properties_name[n], strlen(name)+1, char);
722
strcpy(vars->properties_name[n++], name);
724
strcpy(name, "absorption_coefficient_");
730
if (n != vars->nprop)
731
bft_error(__FILE__, __LINE__, 0,
732
_("number of properties is not correct: %i instead of: %i\n"),
736
bft_printf("==>UIRAPR\n");
737
bft_printf("-->nombre de proprietes = %i\n", vars->nprop);
738
for (i=0 ; i<vars->nprop ; i++)
739
bft_printf("-->properties_ipp[%i]: %i propce[%i]: %i "
740
"properties_name[%i]: %s\n",
741
i, vars->properties_ipp[i],
743
i, vars->properties_name[i]);
747
/*----------------------------------------------------------------------------
748
* Copy variable name from Fortran to C
749
*----------------------------------------------------------------------------*/
751
void CS_PROCF(fcnmra, FCNMRA)
753
const char *const fstr, /* --> Fortran string */
754
int *const len, /* --> String Length */
755
int *const var_id /* --> Variable Id (1 to n) */
764
/* Resize array if necessary */
766
if (*var_id > _cs_gui_max_vars) {
768
if (_cs_gui_max_vars == 0)
769
_cs_gui_max_vars = 16;
771
while (_cs_gui_max_vars <= *var_id)
772
_cs_gui_max_vars *= 2;
774
BFT_REALLOC(_cs_gui_var_rayt, _cs_gui_max_vars, char *);
775
for (i = _cs_gui_last_var; i < _cs_gui_max_vars; i++)
776
_cs_gui_var_rayt[i] = NULL;
779
/* Compute string length (removing start or end blanks) */
782
i1 < *len && (fstr[i1] == ' ' || fstr[i1] == '\t');
786
i2 > i1 && (fstr[i2] == ' ' || fstr[i2] == '\t');
791
/* Should be called once per variable only */
792
assert(_cs_gui_var_rayt[*var_id - 1] == NULL);
796
/* Allocate and copy */
797
BFT_MALLOC(cstr, l + 1, char);
799
for (i = 0 ; i < l ; i++, i1++)
804
_cs_gui_var_rayt[*var_id - 1] = cstr;
808
/* Update variable counter */
809
_cs_gui_last_var = *var_id;
813
/*----------------------------------------------------------------------------
814
* Copy variable name from C to Fortran
815
*----------------------------------------------------------------------------*/
817
void CS_PROCF(cfnmra, CFNMRA)
819
char *const fstr, /* --> Fortran string */
820
int *const len, /* --> String Length */
821
int *const var_id /* --> Variable Id (1 to n) */
829
/* Check that variable name was set */
831
if (*var_id < 1 || *var_id > _cs_gui_last_var)
832
bft_error(__FILE__, __LINE__, 0,
833
_("Name of variable %i was never set.\n"), *var_id);
837
cstr = _cs_gui_var_rayt[*var_id - 1];
841
/* Compute string length (removing start or end blanks) */
847
for (i = 0; i < l; i++)
852
/* Pad with blanks if necessary */
854
for (i = l; i < *len; i++)
858
/*----------------------------------------------------------------------------
859
* Radiative transfer model usray2.F
860
*----------------------------------------------------------------------------*/
862
void CS_PROCF (uiray2, UIRAY2)
864
const int *const itypfb,
865
const int *const iparoi,
866
const int *const iparug,
867
const int *const ivart,
870
const int *const itpimp,
871
const int *const ipgrno,
872
const int *const iprefl,
873
const int *const ifgrno,
874
const int *const ifrefl,
875
const int *const nozppm,
876
const int *const nfabor,
877
const int *const nvar,
887
int output_zone_max = 0;
892
int *faces_list = NULL;
899
zones = cs_gui_boundary_zones_number();
900
output_zone_max = _radiative_boundary_output_zone_max();
902
/* Fisrt iteration only : memory allocation */
903
if (boundary == NULL) {
905
BFT_MALLOC(boundary, 1, cs_radiative_boundary_t);
906
BFT_MALLOC(boundary->label, zones, char* );
907
BFT_MALLOC(boundary->nature, zones, char* );
908
BFT_MALLOC(boundary->output_zone, zones, int );
909
BFT_MALLOC(boundary->type, zones, int );
910
BFT_MALLOC(boundary->emissivity, zones, double );
911
BFT_MALLOC(boundary->thickness, zones, double );
912
BFT_MALLOC(boundary->thermal_conductivity, zones, double );
913
BFT_MALLOC(boundary->external_temp, zones, double );
914
BFT_MALLOC(boundary->internal_temp, zones, double );
915
BFT_MALLOC(boundary->conduction_flux, zones, double );
917
for (izone = 0; izone < zones; izone++) {
919
/* nature, label and description (color or group)
920
of the ith initialization zone */
922
ith_zone = izone + 1;
924
nature = cs_gui_boundary_zone_nature(ith_zone);
926
label = cs_gui_boundary_zone_label(ith_zone);
928
BFT_MALLOC(boundary->label[izone], strlen(label)+1, char);
929
strcpy(boundary->label[izone], label);
931
BFT_MALLOC(boundary->nature[izone], strlen(nature)+1, char);
932
strcpy(boundary->nature[izone], nature);
934
faces_list = cs_gui_get_faces_list(izone,
935
boundaries->label[izone],
936
*nfabor, *nozppm, &faces);
938
/* Default initialization: these values are the same that in raycli
939
but given on each face in raycli whereas here one does not
940
necessarily have boundary faces (parallism) -> duplication */
941
boundary->type[izone] = -1;
942
boundary->output_zone[izone] = -1;
943
boundary->emissivity[izone] = -1.e12;
944
boundary->thickness[izone] = -1.e12;
945
boundary->thermal_conductivity[izone] = -1.e12;
946
boundary->external_temp[izone] = -1.e12;
947
boundary->internal_temp[izone] = -1.e12;
948
boundary->conduction_flux[izone] = 1.e30;
950
if (cs_gui_strcmp(nature, "wall")) {
951
boundary->type[izone] = _radiative_boundary_type(label,
952
*itpimp, *ipgrno, *iprefl,
954
tmp = (double) boundary->output_zone[izone];
955
_radiative_boundary(label, "output_zone", &tmp);
956
boundary->output_zone[izone] = (int) tmp;
957
_radiative_boundary(label, "emissivity", &boundary->emissivity[izone]);
958
_radiative_boundary(label, "thickness", &boundary->thickness[izone]);
959
_radiative_boundary(label, "thermal_conductivity", &boundary->thermal_conductivity[izone]);
960
_radiative_boundary(label, "external_temperature_profile", &boundary->external_temp[izone]);
961
_radiative_boundary(label, "internal_temperature_profile", &boundary->internal_temp[izone]);
962
_radiative_boundary(label, "flux", &boundary->conduction_flux[izone]);
964
} /* if (cs_gui_strcmp(nature, "wall")) */
971
} /* if (boundaries == NULL)*/
973
for (izone = 0; izone < zones; izone++) {
975
/* list of faces building */
978
description = cs_gui_boundary_zone_localization(boundary->label[izone]);
980
fvm_selector_get_list(cs_glob_mesh->select_b_faces,
985
BFT_FREE(description);
988
faces_list = cs_gui_get_faces_list(izone,
989
boundaries->label[izone],
990
*nfabor, *nozppm, &faces);
992
if (cs_gui_strcmp(boundary->nature[izone], "wall"))
994
for (n = 0; n < faces; n++)
996
ifbr = faces_list[n]-1;
998
if (itypfb[ifbr] != *iparoi && itypfb[ifbr] != *iparug)
999
bft_error(__FILE__, __LINE__, 0,
1000
_("One tries to define radiative boundary conditions on boundary which is not a wall.\n"
1001
"The definition of the boundaries natures given in GUI (wall, inlet, outlet,...) \n"
1002
"is modified in a users subroutine (like USCLIM, USCPCL,...). \n"
1003
"The radiative boundary conditions given in GUI must be coherent \n"
1004
"with these new natures.\n"));
1006
izfrdp[ifbr] = boundary->output_zone[izone];
1007
isothp[ifbr] = boundary->type[izone];
1008
if (isothp[ifbr] == *itpimp)
1010
epsp[ifbr] = boundary->emissivity[izone];
1011
tintp[ifbr] = boundary->internal_temp[izone];
1013
else if (isothp[ifbr] == *ipgrno)
1015
xlamp[ifbr] = boundary->thermal_conductivity[izone];
1016
epap[ifbr] = boundary->thickness[izone];
1017
textp[ifbr] = boundary->external_temp[izone];
1018
tintp[ifbr] = boundary->internal_temp[izone];
1019
epsp[ifbr] = boundary->emissivity[izone];
1020
if (boundary->emissivity[izone] == 0.)
1021
isothp[ifbr] = *iprefl;
1023
else if (isothp[ifbr] == *ifgrno)
1025
rcodcl[2 * (*nfabor) * (*nvar) + (*ivart) * (*nfabor) + ifbr] = boundary->conduction_flux[izone];
1026
tintp[ifbr] = boundary->internal_temp[izone];
1027
epsp[ifbr] = boundary->emissivity[izone];
1028
if (boundary->emissivity[izone] != 0.)
1029
isothp[ifbr] = *ifrefl;
1034
j = output_zone_max++;
1035
for (n = 0; n < faces; n++) {
1036
ifbr = faces_list[n]-1;
1039
} /* if nature == "wall" */
1043
BFT_FREE(faces_list);
1046
for (n = 0; n < *nfabor; n++) {
1047
if (izfrdp[n] == -1) iok = 1;
1050
bft_printf("Warning: radiative boundary conditions in GUI are not totally defined \n");
1052
bft_printf("These are radiative boundary conditions defined in GUI: \n");
1053
for (izone = 0; izone < zones; izone++) {
1054
bft_printf(" nature: %s label: %s\n", boundary->nature[izone], boundary->label[izone]);
1055
if (cs_gui_strcmp(boundary->nature[izone], "wall")) {
1056
bft_printf(" output_zone = %i\n", boundary->output_zone[izone]);
1057
bft_printf(" type = %i\n", boundary->type[izone]);
1058
bft_printf(" emissivity = %f\n", boundary->emissivity[izone]);
1059
bft_printf(" thickness= %f\n", boundary->thickness[izone]);
1060
bft_printf(" thermal_conductivity = %f\n", boundary->thermal_conductivity[izone]);
1061
bft_printf(" external_temp = %f\n", boundary->external_temp[izone]);
1062
bft_printf(" internal_temp = %f\n", boundary->internal_temp[izone]);
1063
bft_printf(" conduction_flux= %f\n", boundary->conduction_flux[izone]);
1069
bft_printf("==>UIRAY2\n");
1070
for (izone = 0; izone < zones; izone++) {
1071
bft_printf("--label zone = %s\n", boundary->label[izone]);
1072
if (cs_gui_strcmp(boundary->nature[izone], "wall")) {
1073
bft_printf("----output_zone = %i\n", boundary->output_zone[izone]);
1074
bft_printf("----type = %i\n", boundary->type[izone]);
1075
bft_printf("----emissivity = %f\n", boundary->emissivity[izone]);
1076
bft_printf("----thickness= %f\n", boundary->thickness[izone]);
1077
bft_printf("----thermal_conductivity = %f\n", boundary->thermal_conductivity[izone]);
1078
bft_printf("----external_temp = %f\n", boundary->external_temp[izone]);
1079
bft_printf("----internal_temp = %f\n", boundary->internal_temp[izone]);
1080
bft_printf("----conduction_flux= %f\n", boundary->conduction_flux[izone]);
1086
/*----------------------------------------------------------------------------
1087
* Radiative transfer model usray3.F
1088
*----------------------------------------------------------------------------*/
1091
void CS_PROCF (uiray3, UIRAY3) ( double *const ck,
1092
const int *const ncel,
1098
if (!cs_gui_get_activ_thermophysical_model())
1100
_radiative_transfer_type("absorption_coefficient", &type);
1101
_radiative_transfer_double("absorption_coefficient", &value);
1105
for(i = 0; i < *ncel; i++)
1113
bft_printf("==>UIRAY3\n");
1114
bft_printf("--absorption coefficient type: %d\n", type);
1115
bft_printf("--absorption coefficient by modak: %i\n", imodak);
1117
bft_printf("--absorption coefficient value = %f\n", value);
1122
/*-----------------------------------------------------------------------------
1123
* Free memory: clean global private variables.
1125
* Fortran Interface:
1130
*----------------------------------------------------------------------------*/
1132
void CS_PROCF (memui2, MEMUI2) (void)
1137
if (boundary != NULL) {
1139
/* clean memory for global private structure boundaries */
1141
zones = cs_gui_boundary_zones_number();
1142
for (i=0 ; i < zones ; i++) {
1143
BFT_FREE(boundary->label[i]);
1144
BFT_FREE(boundary->nature[i]);
1146
BFT_FREE(boundary->label);
1147
BFT_FREE(boundary->nature);
1148
BFT_FREE(boundary->output_zone);
1149
BFT_FREE(boundary->type);
1150
BFT_FREE(boundary->emissivity);
1151
BFT_FREE(boundary->thickness);
1152
BFT_FREE(boundary->thermal_conductivity);
1153
BFT_FREE(boundary->external_temp);
1154
BFT_FREE(boundary->internal_temp);
1155
BFT_FREE(boundary->conduction_flux);
1159
/* clean memory for fortran name of variables */
1161
for (i = 0; i < _cs_gui_max_vars; i++)
1162
BFT_FREE(_cs_gui_var_rayt[i]);
1163
BFT_FREE(_cs_gui_var_rayt);
1167
/*----------------------------------------------------------------------------*/