1
/*============================================================================
2
* Management of the GUI parameters file: mobile mesh
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>
50
#include <bft_timer.h>
52
/*----------------------------------------------------------------------------
54
*----------------------------------------------------------------------------*/
56
#include "fvm_selector.h"
58
/*----------------------------------------------------------------------------
60
*----------------------------------------------------------------------------*/
62
#include "mei_evaluate.h"
64
/*----------------------------------------------------------------------------
66
*----------------------------------------------------------------------------*/
69
#include "cs_gui_util.h"
70
#include "cs_gui_variables.h"
71
#include "cs_gui_boundary_conditions.h"
73
#include "cs_prototypes.h"
75
/*----------------------------------------------------------------------------
76
* Header for the current file
77
*----------------------------------------------------------------------------*/
79
#include "cs_gui_mobile_mesh.h"
81
/*----------------------------------------------------------------------------*/
85
/*=============================================================================
86
* Local Macro Definitions
87
*============================================================================*/
89
/* debugging switch */
92
/*============================================================================
94
*============================================================================*/
96
/*----------------------------------------------------------------------------
97
* ALE property choice possible value
98
*----------------------------------------------------------------------------*/
100
enum ale_property_choice
102
ale_property_choice_user_function,
103
ale_property_choice_user_subroutine
106
/*-----------------------------------------------------------------------------
107
* Possible values for boundary nature
108
*----------------------------------------------------------------------------*/
110
enum ale_boundary_nature
112
ale_boundary_nature_fixed_wall,
113
ale_boundary_nature_sliding_wall,
114
ale_boundary_nature_internal_coupling,
115
ale_boundary_nature_external_coupling,
116
ale_boundary_nature_fixed_velocity,
117
ale_boundary_nature_fixed_displacement
120
/*============================================================================
121
* Private function definitions
122
*============================================================================*/
124
/*-----------------------------------------------------------------------------
125
* Return value for iale method
128
* param --> iale parameter
129
* keyword <-- value of the iale parameter
130
*----------------------------------------------------------------------------*/
133
cs_gui_iale_parameter(const char *const param,
134
double *const keyword)
140
path = cs_xpath_init_path();
141
cs_xpath_add_elements(&path, 3, "thermophysical_models", "ale_method", param);
143
if (cs_gui_strcmp(param,"mesh_viscosity")) {
145
cs_xpath_add_attribute(&path, "type");
146
type = cs_gui_get_attribute_value(path);
147
if(cs_gui_strcmp(type, "isotrop"))
149
else if (cs_gui_strcmp(type, "orthotrop"))
152
bft_error(__FILE__, __LINE__, 0, _("Invalid xpath: %s\n"), path);
156
cs_xpath_add_function_text(&path);
157
if (cs_gui_get_double(path, &result)) *keyword = result;
164
/*-----------------------------------------------------------------------------
165
* Return the status of ALE method
168
* keyword <-- status of ale balise
169
*----------------------------------------------------------------------------*/
172
cs_gui_get_ale_status(int *const keyword)
177
path = cs_xpath_init_path();
178
cs_xpath_add_elements(&path, 2, "thermophysical_models", "ale_method");
179
cs_xpath_add_attribute(&path, "status");
181
if(cs_gui_get_status(path, &result))
190
/*-----------------------------------------------------------------------------
191
* Return the viscosity's type of ALE method
194
* type <-- type of viscosity's type
195
*----------------------------------------------------------------------------*/
198
cs_gui_get_ale_viscosity_type(int * type)
203
path = cs_xpath_init_path();
204
cs_xpath_add_elements(&path, 3, "thermophysical_models", "ale_method", "mesh_viscosity");
205
cs_xpath_add_attribute(&path, "type");
207
buff = cs_gui_get_attribute_value(path);
209
if (cs_gui_strcmp(buff, "orthotrop"))
211
else if (cs_gui_strcmp(buff, "isotrop"))
214
bft_error(__FILE__, __LINE__, 0, _("Invalid xpath: %s\n"), path);
221
/*============================================================================
222
* Public function definitions
223
*============================================================================*/
225
/*-----------------------------------------------------------------------------
226
* Initialize mei tree and check for symbols existence
229
* formula --> mei formula
230
* symbols --> array of symbol to check
231
* symbol_nbr --> number of symbol in symbols
232
* variables --> variables required in the formula
233
* variable_nbr --> number of variable in variables
234
* dtref --> time step
235
* ttcabs --> current time
236
* ntcabs --> current iteration number
237
*----------------------------------------------------------------------------*/
240
cs_gui_init_mei_tree(char *formula,
241
const char **symbols,
242
unsigned int symbol_nbr,
243
const char **variables,
244
const double *variables_value,
245
unsigned int variable_nbr,
252
/* return an empty interpreter */
253
mei_tree_t *tree = mei_tree_new(formula);
255
/* Insert variables into mei_tree */
256
for (i = 0; i < variable_nbr; ++i)
260
/* Read value from variables_value if it is not null 0 otherwise */
262
value = variables_value[i];
263
mei_tree_insert(tree, variables[i], value);
266
/* Add commun variables: dt, t, nbIter */
267
mei_tree_insert(tree, "dt", dtref);
268
mei_tree_insert(tree, "t", ttcabs);
269
mei_tree_insert(tree, "iter", ntcabs);
271
/* try to build the interpreter */
272
if (mei_tree_builder(tree))
273
bft_error(__FILE__, __LINE__, 0,
274
_("Error: can not interpret expression: %s\n"), tree->string);
276
/* Check for symbols */
277
for (i = 0; i < symbol_nbr; ++i)
279
const char* symbol = symbols[i];
281
if (mei_tree_find_symbol(tree, symbol))
283
bft_error(__FILE__, __LINE__, 0,
284
_("Error: can not find the required symbol: %s\n"), symbol);
291
/*----------------------------------------------------------------------------
292
* Get the ale property choice
293
*----------------------------------------------------------------------------*/
295
static enum ale_property_choice
296
get_ale_property_choice(void)
300
enum ale_property_choice choice = ale_property_choice_user_function;
301
char *path = cs_xpath_init_path();
303
cs_xpath_add_elements(&path, 3, "thermophysical_models", "ale_method", "property");
304
cs_xpath_add_test_attribute(&path, "label", "mesh_vi1");
305
cs_xpath_add_attribute(&path, "choice");
307
choice_str = cs_gui_get_attribute_value(path);
309
if (cs_gui_strcmp(choice_str , "user_function"))
310
choice = ale_property_choice_user_function;
311
else if (cs_gui_strcmp(choice_str , "user_subroutine"))
312
choice = ale_property_choice_user_subroutine;
314
bft_error(__FILE__, __LINE__, 0,
315
_("Unknow ale property choice %s.\n"), choice);
316
BFT_FREE(choice_str);
321
/*-----------------------------------------------------------------------------
322
* Return the ale boundary nature
323
*----------------------------------------------------------------------------*/
326
get_ale_formula(void)
330
char *path = cs_xpath_short_path();
331
cs_xpath_add_element(&path, "ale_method");
332
cs_xpath_add_element(&path, "formula");
333
cs_xpath_add_function_text(&path);
335
aleFormula = cs_gui_get_text_value(path);
340
/*-----------------------------------------------------------------------------
341
* Return the ale mesh viscosity
342
*----------------------------------------------------------------------------*/
345
get_ale_mesh_viscosity(void)
349
char *path = cs_xpath_short_path();
350
cs_xpath_add_element(&path, "ale_method");
351
cs_xpath_add_element(&path, "mesh_viscosity");
352
cs_xpath_add_attribute(&path, "type");
354
viscosityType = cs_gui_get_attribute_value(path);
356
return viscosityType;
359
/*-----------------------------------------------------------------------------
360
* Get the ale boundary formula
363
* label --> boundary label
364
* choice --> nature: "fixed_velocity" or "fixed_displacement"
365
*----------------------------------------------------------------------------*/
368
get_ale_boundary_formula(const char *const label,
369
const char *const choice)
373
char *path = cs_xpath_init_path();
374
cs_xpath_add_elements(&path, 2, "boundary_conditions", "wall");
375
cs_xpath_add_test_attribute(&path, "label", label);
376
cs_xpath_add_element(&path, "ale");
377
cs_xpath_add_test_attribute(&path, "choice", choice);
378
cs_xpath_add_element(&path, "formula");
379
cs_xpath_add_function_text(&path);
381
formula = cs_gui_get_text_value(path);
387
/*-----------------------------------------------------------------------------
388
* Get uialcl data for fixed displacement
391
* label --> boundary label
392
* begin --> begin index for nodfbr
393
* end --> end index for nodfbr
394
* nnod --> number of node
398
* dtref --> time step
399
* ttcabs --> current time
400
* ntcabs --> current iteration number
401
*----------------------------------------------------------------------------*/
404
uialcl_fixed_displacement(const char *const label,
407
const int *const nnod,
408
const int *const nodfbr,
410
double *const depale,
417
double X_mesh, Y_mesh, Z_mesh;
419
const char* variables[3] = {"mesh_x", "mesh_y", "mesh_z"};
420
unsigned int variable_nbr = 3;
423
char* formula = get_ale_boundary_formula(label, "fixed_displacement");
426
bft_error(__FILE__, __LINE__, 0,
427
_("Boundary nature formula is null for %s.\n"), label);
430
ev = cs_gui_init_mei_tree(formula, variables, variable_nbr,
431
0, 0, 0, dtref, ttcabs, ntcabs);
435
/* Get mei results */
436
X_mesh = mei_tree_lookup(ev, "mesh_x");
437
Y_mesh = mei_tree_lookup(ev, "mesh_y");
438
Z_mesh = mei_tree_lookup(ev, "mesh_z");
441
mei_tree_destroy(ev);
443
/* Set depale and impale */
444
for (ii = begin; ii < end; ++ii)
446
int inod = nodfbr[ii-1] - 1;
447
if (impale[inod] == 0)
449
depale[inod + 0 * (*nnod)] = X_mesh;
450
depale[inod + 1 * (*nnod)] = Y_mesh;
451
depale[inod + 2 * (*nnod)] = Z_mesh;
457
/*-----------------------------------------------------------------------------
458
* Get uialcl data for fixed velocity
461
* label --> boundary label
465
* nfabor --> Number of boundary faces
468
* dtref --> time step
469
* ttcabs --> current time
470
* ntcabs --> current iteration number
471
*----------------------------------------------------------------------------*/
474
uialcl_fixed_velocity(const char* label,
480
double *const rcodcl,
486
const char* variables[3] = { "mesh_u", "mesh_v", "mesh_w" };
489
char* formula = get_ale_boundary_formula(label, "fixed_velocity");
492
bft_error(__FILE__, __LINE__, 0,
493
_("Boundary nature formula is null for %s.\n"), label);
496
ev = cs_gui_init_mei_tree(formula, variables, 3, 0, 0, 0,
497
dtref, ttcabs, ntcabs);
502
rcodcl[ (iuma-1) * nfabor + ifbr ] = mei_tree_lookup(ev, "mesh_u");
503
rcodcl[ (ivma-1) * nfabor + ifbr ] = mei_tree_lookup(ev, "mesh_v");
504
rcodcl[ (iwma-1) * nfabor + ifbr ] = mei_tree_lookup(ev, "mesh_w");
507
mei_tree_destroy(ev);
510
/*-----------------------------------------------------------------------------
511
* Return the ale boundary nature
514
* label --> label of boundary zone
515
*----------------------------------------------------------------------------*/
517
static enum ale_boundary_nature
518
get_ale_boundary_nature(const char *const label)
520
char *ale_boundary_nature;
522
enum ale_boundary_nature nature = ale_boundary_nature_fixed_wall;
524
char *path = cs_xpath_init_path();
525
cs_xpath_add_elements(&path, 2, "boundary_conditions", "wall");
527
cs_xpath_add_test_attribute(&path, "label", label);
528
cs_xpath_add_element(&path, "ale");
529
cs_xpath_add_attribute(&path, "choice");
530
ale_boundary_nature = cs_gui_get_attribute_value(path);
532
if (cs_gui_strcmp(ale_boundary_nature, "fixed_boundary"))
533
nature = ale_boundary_nature_fixed_wall;
534
if (cs_gui_strcmp(ale_boundary_nature, "sliding_boundary"))
535
nature = ale_boundary_nature_sliding_wall;
536
else if (cs_gui_strcmp(ale_boundary_nature, "internal_coupling"))
537
nature = ale_boundary_nature_internal_coupling;
538
else if (cs_gui_strcmp(ale_boundary_nature, "external_coupling"))
539
nature = ale_boundary_nature_external_coupling;
540
else if (cs_gui_strcmp(ale_boundary_nature, "fixed_velocity"))
541
nature = ale_boundary_nature_fixed_velocity;
542
else if (cs_gui_strcmp(ale_boundary_nature, "fixed_displacement"))
543
nature = ale_boundary_nature_fixed_displacement;
546
BFT_FREE(ale_boundary_nature);
551
/*-----------------------------------------------------------------------------
552
* Get boundary attribute like nature or label
555
* ith_zone --> boundary index
556
* nodeName --> xml attribute name. for example, "nature" or "label"
557
*----------------------------------------------------------------------------*/
560
get_boundary_attribute(unsigned int ith_zone,
561
const char *nodeName)
564
char *path = cs_xpath_init_path();
566
cs_xpath_add_element(&path, "boundary_conditions");
567
cs_xpath_add_element_num(&path, "boundary", ith_zone);
568
cs_xpath_add_attribute(&path, nodeName);
570
result = cs_gui_get_attribute_value(path);
577
/*-----------------------------------------------------------------------------
578
* Init xpath for internal coupling with:
580
* boundary_conditions/wall[label=label]/
581
* ale[choice=internal_coupling]/node_name, node_sub_name/text()
584
* label --> boundary label
585
* node_name --> xml node name ("initial_displacement")
586
* node_sub_name --> xml child node of node_name ("X")
587
*----------------------------------------------------------------------------*/
590
init_internal_coupling_xpath(const char* label,
591
const char* node_name,
592
const char* node_sub_name)
596
path = cs_xpath_init_path();
597
cs_xpath_add_elements(&path, 2, "boundary_conditions", "wall");
598
cs_xpath_add_test_attribute(&path, "label", label);
599
cs_xpath_add_element(&path, "ale");
600
cs_xpath_add_test_attribute(&path, "choice", "internal_coupling");
601
cs_xpath_add_element(&path, node_name);
602
cs_xpath_add_element(&path, node_sub_name);
603
cs_xpath_add_function_text(&path);
609
/*-----------------------------------------------------------------------------
610
* Get internal coupling double
613
* label --> boundary label
614
* node_name --> xml node name ("initial_displacement")
615
* node_sub_name --> xml child node of node_name ("X")
616
*----------------------------------------------------------------------------*/
619
get_internal_coupling_double(const char* label,
620
const char* node_name,
621
const char* node_sub_name)
624
char *path = init_internal_coupling_xpath(label, node_name, node_sub_name);
626
if (!cs_gui_get_double(path, &value))
628
bft_error(__FILE__, __LINE__, 0,
629
_("cannot get value for %s %s %s"),
630
label, node_name, node_sub_name);
638
/*-----------------------------------------------------------------------------
639
* Get internal coupling string
642
* label --> boundary label
643
* node_name --> xml node name ("initial_displacement")
644
* node_sub_name --> xml child node of node_name ("formula")
645
*----------------------------------------------------------------------------*/
648
get_internal_coupling_string(const char* label,
649
const char* node_name,
650
const char* node_sub_name)
652
char *path = init_internal_coupling_xpath(label, node_name, node_sub_name);
653
char* str = cs_gui_get_text_value(path);
661
/*-----------------------------------------------------------------------------
662
* Retreive internal coupling x, y and z XML values
665
* label --> boundary label
666
* node_name --> xml node name ("initial_displacement")
667
* xyz <-- result matrix
668
*----------------------------------------------------------------------------*/
671
get_internal_coupling_xyz_values(const char *label,
672
const char *node_name,
675
xyz[0] = get_internal_coupling_double(label, node_name, "X");
676
xyz[1] = get_internal_coupling_double(label, node_name, "Y");
677
xyz[2] = get_internal_coupling_double(label, node_name, "Z");
681
/*-----------------------------------------------------------------------------
682
* Retreive internal coupling advanced windows double value
685
* node_name --> xml node name ("displacement_prediction_alpha")
686
*----------------------------------------------------------------------------*/
689
get_uistr1_advanced_double(const char *const keyword,
693
char *path = cs_xpath_init_path();
695
cs_xpath_add_elements(&path, 3, "thermophysical_models", "ale_method", keyword);
696
cs_xpath_add_function_text(&path);
698
if (cs_gui_get_double(path, &result))
703
/*-----------------------------------------------------------------------------
704
* Retreive internal coupling advanced windows checkbox value
707
* node_name --> xml node name ("monitor_point_synchronisation")
708
*----------------------------------------------------------------------------*/
711
get_uistr1_advanced_checkbox(const char *const keyword, int *const value)
714
char *path = cs_xpath_init_path();
716
cs_xpath_add_elements(&path, 3, "thermophysical_models", "ale_method", keyword);
717
cs_xpath_add_attribute(&path, "status");
719
if (cs_gui_get_status(path, &result))
724
/*-----------------------------------------------------------------------------
725
* Retreive data the internal coupling matrices
728
* label --> boundary label
729
* node_name --> xml matrix node name
730
* symbols --> see cs_gui_init_mei_tree
731
* symbol_nbr --> see cs_gui_init_mei_tree
732
* variables --> see cs_gui_init_mei_tree
733
* variables_value --> see cs_gui_init_mei_tree
734
* variable_nbr --> see cs_gui_init_mei_tree
735
* output_matrix, <-- result matrix
736
* dtref --> time step
737
* ttcabs --> current time
738
* ntcabs --> current iteration number
739
*----------------------------------------------------------------------------*/
742
get_internal_coupling_matrix(const char *label,
743
const char *node_name,
744
const char *symbols[],
745
unsigned int symbol_nbr,
746
const char **variables,
747
const double *variables_value,
748
unsigned int variable_nbr,
749
double *output_matrix,
754
/* Get the formula */
758
char *matrix = get_internal_coupling_string(label, node_name, "formula");
761
bft_error(__FILE__, __LINE__, 0,
762
_("Formula is null for %s %s"), label, node_name);
765
tree = cs_gui_init_mei_tree(matrix, symbols, symbol_nbr,
766
variables, variables_value, variable_nbr,
767
dtref, ttcabs, ntcabs);
770
/* Read matrix values */
771
for (i = 0; i < symbol_nbr; ++i)
773
const char *symbol = symbols[i];
774
output_matrix[i] = mei_tree_lookup(tree, symbol);
777
mei_tree_destroy(tree);
780
/*-----------------------------------------------------------------------------
781
* Retreive data for internal coupling for a specific boundary
784
* label --> boundary label
785
* xmstru <-- Mass matrix
786
* xcstr <-- Damping matrix
787
* xkstru <-- Stiffness matrix
788
* forstr <-- Fluid force matrix
789
* istruc --> internal coupling boundary index
790
* dtref --> time step
791
* ttcabs --> current time
792
* ntcabs --> current iteration number
793
*----------------------------------------------------------------------------*/
796
get_uistr2_data(const char *label,
797
double *const xmstru,
798
double *const xcstru,
799
double *const xkstru,
800
double *const forstr,
806
const char *m_symbols[] = {"m11", "m12", "m13",
808
"m31", "m32", "m33"};
809
const char *c_symbols[] = {"c11", "c12", "c13",
811
"c31", "c32", "c33"};
812
const char *k_symbols[] = {"k11", "k12", "k13",
814
"k31", "k32", "k33"};
816
unsigned int symbol_nbr = sizeof(m_symbols) / sizeof(m_symbols[0]);
818
const char *force_symbols[] = {"fx", "fy", "fz"};
819
unsigned int force_symbol_nbr = sizeof(force_symbols) / sizeof(force_symbols[0]);
821
const unsigned int variable_nbr = 3;
822
const char *variables[3] = {"fluid_fx", "fluid_fy", "fluid_fz"};
823
double variable_values[3];
825
/* Get mass matrix, damping matrix and stiffness matrix */
827
get_internal_coupling_matrix(label, "mass_matrix", m_symbols,
829
&xmstru[istruc * symbol_nbr],
830
dtref, ttcabs, ntcabs);
832
get_internal_coupling_matrix(label, "damping_matrix", c_symbols,
834
&xcstru[istruc * symbol_nbr],
835
dtref, ttcabs, ntcabs);
837
get_internal_coupling_matrix(label, "stiffness_matrix", k_symbols,
839
&xkstru[istruc * symbol_nbr],
840
dtref, ttcabs, ntcabs);
842
/* Set variable for fluid force matrix */
843
variable_values[0] = forstr[istruc * force_symbol_nbr + 0];
844
variable_values[1] = forstr[istruc * force_symbol_nbr + 1];
845
variable_values[2] = forstr[istruc * force_symbol_nbr + 2];
847
/* Get fluid force matrix */
848
get_internal_coupling_matrix(label, "fluid_force_matrix",
849
force_symbols, force_symbol_nbr,
850
variables, variable_values, variable_nbr,
851
&forstr[istruc * force_symbol_nbr],
852
dtref, ttcabs, ntcabs);
855
/*-----------------------------------------------------------------------------
856
* Return the external coupling DDL value
858
* <boundary_conditions>
859
* <wall label=label_argument">
860
* <ale choice="external_coupling">
861
* <node_name_argument choice="off"/>
863
* label --> The wall label.
864
* node_name --> Node name: DDLX, DDLY or DDLZ.
865
*----------------------------------------------------------------------------*/
868
get_external_coupling_ddl(const char *const label,
869
const char *const node_name)
874
char *path = cs_xpath_init_path();
875
cs_xpath_add_elements(&path, 2, "boundary_conditions", "wall");
876
cs_xpath_add_test_attribute(&path, "label", label);
877
cs_xpath_add_element(&path, "ale");
878
cs_xpath_add_test_attribute(&path, "choice", "external_coupling");
879
cs_xpath_add_element(&path, node_name);
880
cs_xpath_add_attribute(&path, "choice");
882
choice = cs_gui_get_attribute_value(path);
883
isOn = cs_gui_strcmp(choice, "on");
891
/*============================================================================
892
* Public Fortran function definitions
893
*============================================================================*/
895
/*----------------------------------------------------------------------------
898
* ncel --> number of cells whithout halo
902
* xyzcen --> cell's gravity center
903
* dtref --> time step
904
* ttcabs --> current time
905
* ntcabs --> current iteration number
906
*----------------------------------------------------------------------------*/
908
void CS_PROCF (uivima, UIVIMA) (const cs_int_t *const ncel,
909
double *const viscmx,
910
double *const viscmy,
911
double *const viscmz,
912
const double *const xyzcen,
914
double *const ttcabs,
915
const int *const ntcabs)
918
const char* symbols[3] = { "x", "y", "z" };
919
const char* variables[3] = { "mesh_vi1", "mesh_vi2", "mesh_vi3" };
920
unsigned int variable_nbr = 1;
922
/* Check if user function is selection */
923
if (get_ale_property_choice() == ale_property_choice_user_function)
927
char *aleFormula = get_ale_formula();
928
char *viscosityType = get_ale_mesh_viscosity();
929
unsigned int isOrthotrop = cs_gui_strcmp(viscosityType, "orthotrop");
935
bft_error(__FILE__, __LINE__, 0,
936
_("Formula is null for ale.\n"));
939
ev = cs_gui_init_mei_tree(aleFormula, variables,
940
variable_nbr, symbols, 0, 3,
941
*dtref, *ttcabs, *ntcabs);
943
/* for each cell, update the value of the table of symbols for each scalar
944
(including the thermal scalar), and evaluate the interpreter */
945
for (iel = 0; iel < *ncel; iel++)
948
mei_tree_insert(ev, "x", xyzcen[3 * iel + 0]);
949
mei_tree_insert(ev, "y", xyzcen[3 * iel + 1]);
950
mei_tree_insert(ev, "z", xyzcen[3 * iel + 2]);
954
/* Set viscmx, viscmy and viscmz */
955
viscmx[iel] = mei_tree_lookup(ev, "mesh_vi1");
958
viscmy[iel] = mei_tree_lookup(ev, "mesh_vi2");
959
viscmz[iel] = mei_tree_lookup(ev, "mesh_vi3");
962
mei_tree_destroy(ev);
963
BFT_FREE(aleFormula);
964
BFT_FREE(viscosityType);
968
/*----------------------------------------------------------------------------
969
* ALE related keywords
976
* INTEGER IALE <-- iale method activation
977
* INTEGER NALINF <-- number of sub iteration of initialization
979
* INTEGER NALIMX <-- max number of iterations of implicitation of
980
* the displacement of the structures
981
* DOUBLE EPALIM <-- realtive precision of implicitation of
982
* the displacement of the structures
983
* INTEGER IORTVM <-- type of viscosity of mesh
985
*----------------------------------------------------------------------------*/
987
void CS_PROCF (uialin, UIALIN) (int *const iale,
990
double *const epalim,
995
cs_gui_get_ale_status(iale);
998
value =(double) *nalinf;
999
cs_gui_iale_parameter("fluid_initialization_sub_iterations", &value);
1000
*nalinf = (int) value;
1002
value =(double) *nalimx;
1003
cs_gui_iale_parameter("max_iterations_implicitation", &value);
1004
*nalimx = (int) value;
1006
cs_gui_iale_parameter("implicitation_precision", epalim);
1008
value =(double) *iortvm;
1009
cs_gui_iale_parameter("mesh_viscosity", &value);
1010
*iortvm = (int) value;
1014
bft_printf("==>UIALIN\n");
1015
bft_printf("--iale = %i\n", *iale);
1017
bft_printf("--nalinf = %i\n", *nalinf);
1018
bft_printf("--nalimx = %i\n", *nalimx);
1019
bft_printf("--epalim = %g\n", *epalim);
1020
bft_printf("--iortvm = %i\n", *iortvm);
1025
/*-----------------------------------------------------------------------------
1029
* nfabor --> Number of boundary faces
1030
* nozppm --> Max number of boundary conditions zone
1032
* ipnfbr --> First node position for each boundary in nodfbr
1033
* nodfbr --> uialcl_fixed_displacement
1034
* impale --> uialcl_fixed_displacement
1035
* depale --> See uialcl_fixed_displacement
1036
* dtref --> time step
1037
* ttcabs --> current time
1038
* ntcabs --> current iteration number
1039
* iuma --> See uialcl_fixed_velocity
1040
* ivma --> See uialcl_fixed_velocity
1041
* iwma --> See uialcl_fixed_velocity
1042
* rcodcl --> See uialcl_fixed_velocity
1043
*----------------------------------------------------------------------------*/
1045
void CS_PROCF (uialcl, UIALCL) (const int *const nfabor,
1046
const int *const nozppm,
1047
const int *const ibfixe,
1048
const int *const igliss,
1049
const int *const ivimpo,
1051
const int *const ipnfbr,
1052
const int *const nnod,
1053
const int *const nodfbr,
1055
double *const depale,
1056
double *const dtref,
1057
double *const ttcabs,
1058
const int *const ntcabs,
1059
const int *const iuma,
1060
const int *const ivma,
1061
const int *const iwma,
1062
double *const rcodcl)
1069
int zones = cs_gui_boundary_zones_number();
1071
/* At each time-step, loop on boundary faces: */
1072
for (izone=0 ; izone < zones ; izone++)
1074
int* faces_list = cs_gui_get_faces_list(izone,
1075
boundaries->label[izone],
1080
/* get the ale choice */
1081
const char* label = boundaries->label[izone];
1082
enum ale_boundary_nature nature = get_ale_boundary_nature(label);
1084
if (nature == ale_boundary_nature_fixed_wall)
1086
for (ifac = 0; ifac < faces; ifac++)
1088
int ifbr = faces_list[ifac]-1;
1089
ialtyb[ifbr] = *ibfixe;
1092
else if (nature == ale_boundary_nature_sliding_wall)
1094
for (ifac = 0; ifac < faces; ifac++)
1096
int ifbr = faces_list[ifac]-1;
1097
ialtyb[ifbr] = *igliss;
1100
else if (nature == ale_boundary_nature_fixed_displacement)
1102
t0 = bft_timer_wtime();
1103
for (ifac = 0; ifac < faces; ifac++)
1105
int ifbr = faces_list[ifac]-1;
1106
uialcl_fixed_displacement(label, ipnfbr[ifbr], ipnfbr[ifbr+1],
1107
nnod, nodfbr, impale, depale,
1108
*dtref, *ttcabs, *ntcabs);
1110
cs_gui_add_mei_time(bft_timer_wtime() - t0);
1112
else if (nature == ale_boundary_nature_fixed_velocity)
1114
t0 = bft_timer_wtime();
1115
for (ifac = 0; ifac < faces; ifac++)
1117
int ifbr = faces_list[ifac]-1;
1118
uialcl_fixed_velocity(label, *iuma, *ivma, *iwma,
1119
*nfabor, ifbr, rcodcl,
1120
*dtref, *ttcabs, *ntcabs);
1121
ialtyb[ifbr] = *ivimpo;
1123
cs_gui_add_mei_time(bft_timer_wtime() - t0);
1125
BFT_FREE(faces_list);
1129
/*-----------------------------------------------------------------------------
1130
* Retreive data for internal coupling. Called once at initialization
1133
* nfabor --> Number of boundary faces
1134
* idfstr --> Structure definition
1135
* aexxst <-- Displacement prediction alpha
1136
* bexxst <-- Displacement prediction beta
1137
* cfopre <-- Stress prediction alpha
1138
* ihistr <-- Monitor point synchronisation
1139
* xstr0 <-- Values of the initial displacement
1140
* xstreq <-- Values of the equilibrium displacement
1141
* vstr0 <-- Values of the initial velocity
1142
*----------------------------------------------------------------------------*/
1144
void CS_PROCF (uistr1, UISTR1) (const int *const nfabor,
1159
int *faces_list = NULL;
1160
unsigned int istruct = 0;
1162
/* Get advanced data */
1163
get_uistr1_advanced_double("displacement_prediction_alpha", aexxst);
1164
get_uistr1_advanced_double("displacement_prediction_beta", bexxst);
1165
get_uistr1_advanced_double("stress_prediction_alpha", cfopre);
1166
get_uistr1_advanced_checkbox("monitor_point_synchronisation", ihistr);
1168
zones = cs_gui_boundary_zones_number();
1170
/* At each time-step, loop on boundary faces */
1171
for (izone=0 ; izone < zones ; izone++)
1173
char *nature = get_boundary_attribute(izone + 1, "nature");
1174
char *label = get_boundary_attribute(izone + 1, "label");
1176
/* Keep only internal coupling */
1177
if (get_ale_boundary_nature(label) == ale_boundary_nature_internal_coupling)
1179
/* Read initial_displacement, equilibrium_displacement and initial_velocity */
1180
get_internal_coupling_xyz_values(label, "initial_displacement",
1181
&xstr0[3 * istruct]);
1182
get_internal_coupling_xyz_values(label, "equilibrium_displacement",
1183
&xstreq[3 * istruct]);
1184
get_internal_coupling_xyz_values(label, "initial_velocity",
1185
&vstr0[3 * istruct]);
1187
faces_list = cs_gui_get_faces_list(izone, label, *nfabor, 0, &faces);
1188
/* Set idfstr to positiv index starting at 1 */
1189
for (ifac = 0; ifac < faces; ifac++)
1191
ifbr = faces_list[ifac]-1;
1192
idfstr[ifbr] = istruct + 1;
1195
BFT_FREE(faces_list);
1202
/*-----------------------------------------------------------------------------
1203
* Retreive data for internal coupling. Called at each step
1206
* xmstru <-- Mass matrix
1207
* xcstr <-- Damping matrix
1208
* xkstru <-- Stiffness matrix
1209
* forstr <-- Fluid force matrix
1210
* dtref --> time step
1211
* ttcabs --> current time
1212
* ntcabs --> current iteration number
1213
*----------------------------------------------------------------------------*/
1215
void CS_PROCF (uistr2, UISTR2) (double *const xmstru,
1216
double *const xcstru,
1217
double *const xkstru,
1218
double *const forstr,
1219
double *const dtref,
1220
double *const ttcabs,
1224
unsigned int istru = 0;
1226
int zones = cs_gui_boundary_zones_number();
1228
/* At each time-step, loop on boundary faces */
1229
for (izone=0 ; izone < zones ; izone++)
1231
const char *label = boundaries->label[izone];
1233
/* Keep only internal coupling */
1234
if (get_ale_boundary_nature(label) == ale_boundary_nature_internal_coupling)
1237
/* Read internal coupling data for boundaries */
1238
for (int ii=0; ii<3; ii++)
1245
get_uistr2_data(label,
1259
/*-----------------------------------------------------------------------------
1260
* Retreive data for external coupling
1263
* nfabor <-- Number of boundary faces
1264
* idfstr <-- Structure definition
1265
* asddlf <-- Block of the DDL forces
1266
*----------------------------------------------------------------------------*/
1269
CS_PROCF (uiaste, UIASTE) (const int *const nfabor,
1271
double *const asddlf)
1279
int zones = cs_gui_boundary_zones_number();
1281
/* At each time-step, loop on boundary faces */
1282
for (izone=0 ; izone < zones ; izone++)
1284
const char *label = boundaries->label[izone];
1286
/* Keep only internal coupling */
1287
if (get_ale_boundary_nature(label) == ale_boundary_nature_external_coupling)
1289
int* faces_list = cs_gui_get_faces_list(izone, label, *nfabor, 0, &faces);
1291
/* Get DDLX, DDLY and DDLZ values */
1292
asddlf[istruct * 3 + 0] = get_external_coupling_ddl(label, "DDLX") ? 0 : 1;
1293
asddlf[istruct * 3 + 1] = get_external_coupling_ddl(label, "DDLY") ? 0 : 1;
1294
asddlf[istruct * 3 + 2] = get_external_coupling_ddl(label, "DDLZ") ? 0 : 1;
1296
/* Set idfstr with negativ value starting from -1 */
1297
for (ifac = 0; ifac < faces; ifac++)
1299
ifbr = faces_list[ifac]-1;
1300
idfstr[ifbr] = -istruct - 1;
1303
BFT_FREE(faces_list);
1308
/*----------------------------------------------------------------------------*/