-
Notifications
You must be signed in to change notification settings - Fork 0
/
bc_colloc.c
2808 lines (2518 loc) · 82.3 KB
/
bc_colloc.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/************************************************************************ *
* Goma - Multiphysics finite element software *
* Sandia National Laboratories *
* *
* Copyright (c) 2014 Sandia Corporation. *
* *
* Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, *
* the U.S. Government retains certain rights in this software. *
* *
* This software is distributed under the GNU General Public License. *
\************************************************************************/
/*
* Routines for calculating Boundary Conditions and adding them
* to matrix_fill
*
*
*$Id: bc_colloc.c,v 5.11 2010-07-21 16:39:26 hkmoffa Exp $
*/
/* Standard include files */
#ifdef HAVE_CONFIG_H
#include <config.h>
#endif
#include <stdio.h>
#include <stdlib.h>
#include <strings.h>
#include <string.h>
#include <math.h>
#include <stdlib.h>
#include "dpi.h"
#define _BC_COLLOC_C
#include "goma.h"
#ifdef USE_CGM
#include "gm_cgm_c_interface.h"
#endif
/*******************************************************************************/
int
apply_point_colloc_bc (
int ija[], /* Vector of integer pointers into the vector a */
double a[], /* Jacobian Matrix */
double resid_vector[], /* Residual vector for the current processor */
double delta_t, /* current time step size */
double theta, /* parameter to vary time integration:
* explicit (theta = 1) -- implicit (theta = 0) */
int ielem, /* element number */
int ip_total, /* total number of gauss points */
int ielem_type, /* element type */
int num_local_nodes,
int ielem_dim,
int iconnect_ptr,
ELEM_SIDE_BC_STRUCT *elem_side_bc,
/* Pointer to an element side boundary
* condition structure */
int num_total_nodes,
int local_node_list_fs[],/* dimensioned [MDE]; list to keep track of
* nodes at which solid contributions have been
* transfered to liquid (fluid-solid
* boundaries) */
double time_value,
Exo_DB *exo)
/*******************************************************************************
Function which fills the FEM stiffness matrices and the right-hand-side vector
with contributions from boundary conditions.
Authors: Harry Moffat, Scott Hutchinson (1421) and others
Transferred to separate files by Rich Cairncross
Date: 12 December 1994
******************************************************************************/
{
int w, i, I, ibc, j, id, err, var, eqn, ldof_eqn, ldof_var, icount;
int ieqn, pvar;
int el1, el2, nf, jk;
int status = 0;
int index_eq, matID_apply;
int bc_input_id;
int var_flag;
int dof_map[MDE];
int n_dof[MAX_VARIABLE_TYPES];
int n_dofptr[MAX_VARIABLE_TYPES][MDE];
int doMeshMapping = 0;
double xi[DIM]; /* Gaussian-quadrature point locations */
double x_dot[MAX_PDIM];
double dsigma_dx[DIM][MDE];
double phi_j;
double func;
double f_time;
double d_func[MAX_VARIABLE_TYPES + MAX_CONC];
double kfunc[DIM];
double d_kfunc[DIM][MAX_VARIABLE_TYPES + MAX_CONC][MDE];
double time_intermediate = time_value-theta*delta_t;
/* time at which bc's are evaluated */
const double penalty = BIG_PENALTY;
VARIABLE_DESCRIPTION_STRUCT *vd;
int doFullJac = 0;
double nwall[3];
int contact_flag = FALSE;
/***************************************************************************/
/* START OF SURFACE LOOPS THAT DON'T REQUIRE INTEGRATION (STRONG SENSE) */
/* (POINTWISE COLLOcATION) */
/* LOOP OVER THE SURFACES IN THE CURRENT ELEMENT THAT */
/* NEED THE APPLICATION OF A STONG SENSE CONDITION */
/* - INITIALIZATION THAT IS DEPENDENT ON THE IDENTITY OF THE SURFACE */
/***************************************************************************/
/***************************************************************************/
/* LOOP OVER THE LOCAL ELEMENT NODE NUMBER */
/* FOR NODES THAT LIE ON THE CURRENT SURFACE */
/* - INITIALIZATION THAT IS DEPENDENT ON THE LOCAL ELEMENT NODE NUMBER */
/***************************************************************************/
#ifdef DEBUG_HKM
if (ei->ielem == 1032) {
// printf("we are here - 1032\n");
}
#endif
for (i = 0; i < (int) elem_side_bc->num_nodes_on_side; i++) {
/* Find the local element node number for the current node */
id = (int) elem_side_bc->local_elem_node_id[i];
/* Find the local node number given the local element node
* number, 'id'
* I is the global node number of this node !!
*/
I = Proc_Elem_Connect[iconnect_ptr + id];
/*
* Check to see if the current local node number is owned by the
* processor - this is done by checking I < num_owned_nodes
*/
if (I < DPI_ptr->num_owned_nodes) {
/*
* Load the field variable values at this node point
*/
find_nodal_stu(id, ielem_type, &xi[0], &xi[1], &xi[2]);
err = load_basis_functions( xi, bfd );
EH( err, "problem from load_basis_functions");
err = beer_belly();
EH( err, "beer_belly");
err = load_fv();
EH( err, "load_fv");
err = load_bf_grad();
EH( err, "load_bf_grad");
err = load_bf_mesh_derivs();
EH( err, "load_bf_mesh_derivs");
/* calculate the shape functions and their gradients */
/* calculate the determinant of the surface jacobian and the normal to
* the surface all at one time */
surface_determinant_and_normal(ielem, iconnect_ptr, num_local_nodes,
ielem_dim - 1,
(int) elem_side_bc->id_side,
(int) elem_side_bc->num_nodes_on_side,
(elem_side_bc->local_elem_node_id) );
if (ielem_dim !=3) {
calc_surf_tangent(ielem, iconnect_ptr, num_local_nodes, ielem_dim-1,
(int) elem_side_bc->num_nodes_on_side,
(elem_side_bc->local_elem_node_id));
}
/*
* Load up porous media variables and properties, if needed
*/
if (mp->PorousMediaType == POROUS_UNSATURATED ||
mp->PorousMediaType == POROUS_TWO_PHASE ||
mp->PorousMediaType == POROUS_SATURATED) {
err = load_porous_properties();
EH(err, "load_porous_properties");
}
if (TimeIntegration != STEADY) {
if (mp->SurfaceTensionModel != CONSTANT) {
load_surface_tension(dsigma_dx);
if (neg_elem_volume) return (status);
}
}
if (TimeIntegration != STEADY) {
for (icount = 0; icount < ielem_dim; icount++ ) {
x_dot[icount] = fv_dot->x[icount];
}
} else {
for (icount = 0; icount < ielem_dim; icount++ ) {
x_dot[icount] = 0.0;
}
}
do_LSA_mods(LSA_SURFACE);
/**********************************************************************/
/* */
/* LOOP OVER THE INDIVIDUAL SURFACE CONDITION TERMS */
/* Calculate the functions for each collocation */
/* condition at this point */
/* */
/**********************************************************************/
for (ibc = 0;
(bc_input_id = (int) elem_side_bc->BC_input_id[ibc]) != -1 ;
ibc++) {
/*
* This function only handles COLLOCATE_SURF boundary conditions.
* All other boundary conditions are weeded out here
*/
if (BC_Types[bc_input_id].desc->method == COLLOCATE_SURF) {
/* Initialize time function to unity */
f_time = 1.0;
/*
* Check to see if this BC on this node is applicable
* (i.e. no other overriding Dirichlet conditions,
* And find the global unknown number for applying this condition
*/
index_eq = bc_eqn_index(id, I, bc_input_id, ei->mn,
0, &eqn, &matID_apply, &vd);
if (index_eq >= 0) {
/* initialize the general function to zero */
func = 0.;
init_vec_value(d_func, 0.0, MAX_VARIABLE_TYPES + MAX_CONC);
doFullJac = 0;
doMeshMapping = 0;
/*
* Here's a RECIPE for adding new boundary conditions so you don't have any
* excuses not to add new ones. The changes should be made in at least
* four files (rf_bc_const.h, mm_names.h, mm_input.c, and bc_[method].c)
* for some boundary conditions you may want to make additions elsewhere also.
* One example of extra additions is in el_exoII_io.c where the ss_dup_list
* is created - you may want to adapt the logic for rotation of new conditions.
* (note that these lines are repeated at each place where you need to
* make changes):
* Boundary Condition Step 1: add Macro Substitution for your new boundary
* condition in rf_bc_const.h - this is how you
* rf_bc_const.h will refer to the new boundary condition
* throughout the code. Make sure the integer
* you choose is unique from all the other BC
* types.
* Boundary Condition Step 2: add a description of your boundary condition
* to the BC_Desc structure array in mm_names.h.
* mm_names.h This structure includes information about the
* type of boundary condition, which equation it
* applies to, what variables it is sensitive to,
* whether to rotate mesh or momentum equations,
* etc. It is very important that you fill out
* this structure carefully, otherwise the code
* won't know what to do.
* Boundary Condition Step 3: add your BC case to the correct format listing
* for reading the necessary arguments from the
* mm_input_bc.c input file in mm_input_bc.c.
*
* Boundary Condition Step 4: Add a function call (and a function) in the
* correct routine for evaluating your boundary
* bc_colloc.c condition. This will probably in bc_colloc.c
* bc_integ.c for collocated conditions or bc_integ.c for
* strong or weak integrated conditions.
* Boundary Condition Step 5: use and enjoy your new boundary condition
*
* Step 4 should be done below (or in bc_integ.c), and add your new function at
* the end of this file (see fplane)
*/
/* load in functional representation of this condition, and its sensitivity
* to all nodal unknowns */
switch (BC_Types[bc_input_id].BC_Name) {
/* Section for surface terms involving the mesh equations */
case DXDISTNG_BC:
case DYDISTNG_BC:
case DZDISTNG_BC:
case DISTNG_BC:
fTmelting (&func, d_func,
BC_Types[bc_input_id].BC_Data_Float[0]);
break;
case SPLINEX_BC:
case SPLINEY_BC:
case SPLINEZ_BC:
case SPLINE_BC:
fspline(ielem_dim, &func, d_func, BC_Types[bc_input_id].u_BC,
time_intermediate);
break;
case SPLINEX_RS_BC:
case SPLINEY_RS_BC:
case SPLINEZ_RS_BC:
case SPLINE_RS_BC:
fspline_rs(ielem_dim, &func, d_func, BC_Types[bc_input_id].u_BC,
time_intermediate);
break;
case T_USER_BC:
tuser(&func, d_func, BC_Types[bc_input_id].u_BC,
time_intermediate);
break;
case DX_USER_BC:
dx_user_surf(&func, d_func, BC_Types[bc_input_id].u_BC, time_intermediate);
break;
case DY_USER_BC:
dy_user_surf(&func, d_func, BC_Types[bc_input_id].u_BC, time_intermediate);
break;
case DZ_USER_BC:
dz_user_surf(&func, d_func, BC_Types[bc_input_id].u_BC, time_intermediate);
break;
case P_LIQ_USER_BC:
p_liq_user_surf(&func, d_func, BC_Types[bc_input_id].u_BC, time_intermediate);
break;
case SH_P_OPEN_USER_BC:
shell_p_open_user_surf(&func, d_func, BC_Types[bc_input_id].u_BC, time_intermediate);
break;
case POROUS_LIQ_PRESSURE_FILL_BC:
{
MATRL_PROP_STRUCT *mp_2=NULL;
porous_liq_fill(&func, d_func, BC_Types[bc_input_id].BC_Data_Int[0],
BC_Types[bc_input_id].BC_Data_Int[1],
BC_Types[bc_input_id].BC_Data_Float[0],
BC_Types[bc_input_id].BC_Data_Float[1],
BC_Types[bc_input_id].BC_Data_Float[2],
mp_2);
}
break;
case PLANEX_BC:
case PLANEY_BC:
case PLANEZ_BC:
case PLANE_BC:
fplane(ielem_dim, &func, d_func,
BC_Types[bc_input_id].BC_Data_Float);
break;
#ifdef USE_CGM
case SM_PLANE_BC: /* Solid Model PLANE BC */
/* I took out the plane generation at the BC level. It
has been delayed. The intention is that it will be
created on it's first call, and then just referred to
thereafter. MMH */
/* sm_fplane (ielem_dim, &func, d_func, */
/* BC_Types[bc_input_id].CGM_plane_handle); */
EH(-1, "CGM::bc_colloc.c:320 not implemented yet");
break;
#endif
case MOVING_PLANE_BC:
{
double t = time_intermediate;
fplane (ielem_dim, &func, d_func,
BC_Types[bc_input_id].BC_Data_Float);
func += (BC_Types[bc_input_id].BC_Data_Float[4]*t +
BC_Types[bc_input_id].BC_Data_Float[5]*t*t +
BC_Types[bc_input_id].BC_Data_Float[6]*t*t*t);
if (af->Assemble_LSA_Mass_Matrix)
EH(-1, "LSA is not currently compatible with MOVING_PLANE_BC");
}
break;
case MESH_CONSTRAINT_BC:
fmesh_constraint(&func, d_func, bc_input_id);
break;
case VVARY_BC:
case UVARY_BC:
case WVARY_BC:
var_flag = BC_Types[bc_input_id].desc->equation;
fvelocity_profile(var_flag, ielem_dim, BC_Types[bc_input_id].BC_Name,
&func, d_func, BC_Types[bc_input_id].u_BC,
time_intermediate);
break;
case GD_CONST_BC:
case GD_LINEAR_BC:
case GD_INVERSE_BC:
case GD_PARAB_BC:
case GD_PARAB_OFFSET_BC:
case GD_CIRC_BC:
case GD_POLYN_BC:
case GD_TABLE_BC:
err = fgeneralized_dirichlet(&func, d_func, BC_Types[bc_input_id].BC_Name,
bc_input_id, theta, delta_t);
EH(err, "Illegal entry in Generalized Dirichlet Condition ");
break;
case GD_TIME_BC:
err = evaluate_time_func(time_intermediate, &f_time, bc_input_id);
EH(err, "Problems in evaluating time function");
if (af->Assemble_LSA_Mass_Matrix)
EH(-1, "LSA is not currently compatible with GD_TIME_BC");
break;
case POROUS_PRESSURE_BC:
porous_pressure(&func, d_func,
(int) BC_Types[bc_input_id].BC_Data_Int[0],
(int) BC_Types[bc_input_id].BC_Data_Int[1]);
break;
case POROUS_PRESSURE_LUB_BC:
porous_pressure_lub(&func, d_func,
elem_side_bc->id_side, xi, exo,
BC_Types[bc_input_id].BC_Data_Float[0]);
break;
case LUBP_SH_FP_FLUX_BC:
put_lub_flux_in_film(id, I, ielem_dim, ija, a, resid_vector,
(int) BC_Types[bc_input_id].BC_Data_Int[0],
(int) BC_Types[bc_input_id].BC_Data_Int[1],
local_node_list_fs);
func = 0.;
break;
case FLUID_SOLID_BC:
case SOLID_FLUID_BC:
put_liquid_stress_in_solid(id, I,
ielem_dim, ija, a , resid_vector,
(int) BC_Types[bc_input_id].BC_Data_Int[0],
(int) BC_Types[bc_input_id].BC_Data_Int[1],
local_node_list_fs,
BC_Types[bc_input_id].BC_Data_Float[0]);
func = 0.; /* this boundary condition rearranges values already in res and jac,
and does not add anything into the residual */
break;
case SOLID_FLUID_RS_BC:
put_liquid_stress_in_solid_ALE(id, I,
ielem_dim, ija, a , resid_vector,
(int) BC_Types[bc_input_id].BC_Data_Int[0],
(int) BC_Types[bc_input_id].BC_Data_Int[1],
local_node_list_fs,
BC_Types[bc_input_id].BC_Data_Float[0]);
func = 0.; /* this boundary condition rearranges values already in res and jac,
* and does not add anything into the residual */
break;
case FLUID_SOLID_RS_BC:
EH(-1,"FLUID_SOLID_RS bc not implemented yet");
break;
case SH_FLUID_STRESS_BC:
/*Note that we send in i, with id, as this is the shell-element counterpart local num */
put_fluid_stress_on_shell(id, i , I,
ielem_dim, ija, a , resid_vector,
local_node_list_fs,
BC_Types[bc_input_id].BC_Data_Float[0]);
func = 0.; /* this boundary condition rearranges values already in res and jac,
* and does not add anything into the residual */
break;
case KINEMATIC_COLLOC_BC:
case VELO_NORM_COLLOC_BC:
/* initialize the general function to zero may have more than
* one entry for vector conditions like capillary */
memset(kfunc, 0, DIM*sizeof(double) );
memset(d_kfunc,0, DIM*(MAX_VARIABLE_TYPES + MAX_CONC)*MDE*sizeof(double));
fvelo_normal_bc(kfunc, d_kfunc,
BC_Types[bc_input_id].BC_Data_Float[0], contact_flag = FALSE,
x_dot, theta, delta_t,
(int) BC_Types[bc_input_id].BC_Name,0,0, 135.0);
doFullJac = 1;
func = kfunc[0];
break;
case VELO_NORMAL_LS_COLLOC_BC:
/* initialize the general function to zero may have more than
* one entry for vector conditions like capillary */
memset(kfunc, 0, DIM*sizeof(double) );
memset(d_kfunc,0, DIM*(MAX_VARIABLE_TYPES + MAX_CONC)*MDE*sizeof(double));
contact_flag = (ls != NULL);
fvelo_normal_bc(kfunc, d_kfunc,
BC_Types[bc_input_id].BC_Data_Float[0], contact_flag,
x_dot, theta, delta_t,
(int) BC_Types[bc_input_id].BC_Name,
BC_Types[bc_input_id].BC_Data_Float[1],
BC_Types[bc_input_id].BC_Data_Float[2],
BC_Types[bc_input_id].BC_Data_Float[3]);
doFullJac = 1;
func = kfunc[0];
break;
case KIN_DISPLACEMENT_COLLOC_BC:
memset(kfunc, 0, DIM*sizeof(double) );
memset(d_kfunc,0, DIM*(MAX_VARIABLE_TYPES + MAX_CONC)*MDE*sizeof(double));
f_kinematic_displacement_bc(kfunc, d_kfunc,
BC_Types[bc_input_id].BC_Data_Int[0],
BC_Types[bc_input_id].BC_ID,
BC_Types[bc_input_id].u_BC,
BC_Types[bc_input_id].len_u_BC);
func = kfunc[0];
doFullJac = 1;
break;
case TABLE_BC:
apply_table_bc(&func, d_func, &BC_Types[bc_input_id], time_value );
if(neg_elem_volume) return(status);
break;
case SH_GAMMA1_DERIV_SYMM_BC:
memset(kfunc, 0, DIM*sizeof(double) );
memset(d_kfunc,0, DIM*(MAX_VARIABLE_TYPES + MAX_CONC)*MDE*sizeof(double));
fgamma1_deriv_bc(kfunc, d_kfunc, BC_Types[bc_input_id].BC_Data_Float[0]);
func = kfunc[0];
doFullJac = 1;
el1 = ei->ielem;
nf = num_elem_friends[el1];
if (nf == 0) {
EH(-1, "no friends");
};
el2 = elem_friends[el1][0];
err = load_neighbor_var_data(el1, el2, n_dof, dof_map, n_dofptr, -1, xi, exo);
doMeshMapping = 1;
break;
case SH_GAMMA2_DERIV_SYMM_BC:
memset(kfunc, 0, DIM*sizeof(double) );
memset(d_kfunc,0, DIM*(MAX_VARIABLE_TYPES + MAX_CONC)*MDE*sizeof(double));
fgamma2_deriv_bc(kfunc, d_kfunc, BC_Types[bc_input_id].BC_Data_Float[0]);
func = kfunc[0];
doFullJac = 1;
el1 = ei->ielem;
nf = num_elem_friends[el1];
if (nf == 0) {
EH(-1, "no friends");
};
el2 = elem_friends[el1][0];
err = load_neighbor_var_data(el1, el2, n_dof, dof_map, n_dofptr, -1, xi, exo);
doMeshMapping = 1;
break;
case DVZDR_ZERO_BC:
memset(kfunc, 0, DIM*sizeof(double) );
memset(d_kfunc,0, DIM*(MAX_VARIABLE_TYPES + MAX_CONC)*MDE*sizeof(double));
nwall[0] = BC_Types[bc_input_id].BC_Data_Float[1];
nwall[1] = BC_Types[bc_input_id].BC_Data_Float[2];
nwall[2] = BC_Types[bc_input_id].BC_Data_Float[3];
dvzdr_zero_deriv_bc(kfunc, d_kfunc, nwall, BC_Types[bc_input_id].BC_Data_Float[0]);
func = kfunc[0];
doFullJac = 1;
break;
} /* end of SWITCH statement */
/************************************************************************/
/* */
/* ADD the Boundary condition functions into the Residual */
/* vector and Jacobian Matrix */
/* */
/************************************************************************/
/*
* Collocated boundary conditions are always applied on the first
* dof at a node. They are not discontinuous variables friendly
*/
ldof_eqn = ei->ln_to_first_dof[eqn][id];
if (eqn == R_MASS) {
ieqn = MAX_PROB_EQN + BC_Types[bc_input_id].species_eq;
} else {
ieqn = upd->ep[eqn];
}
if (ldof_eqn != -1) {
lec->R[ieqn][ldof_eqn] += penalty * func;
lec->R[ieqn][ldof_eqn] *= f_time;
/*
* add sensitivities into matrix
* - find index of sensitivity in matrix
* (if variable is not defined at this node,
* loop over all dofs in element)
* - add into matrix
*/
if (af->Assemble_Jacobian) {
for (var = 0; var < MAX_VARIABLE_TYPES; var++) {
pvar = upd->vp[var];
if (pvar != -1 && (BC_Types[bc_input_id].desc->sens[var] || 1)) {
/*
* Warning!!!!!!!!!!!!!!!!!!!!!!!
* This section of code will need to be revamped
* to work for discontinuous variable.
*/
if (var != MASS_FRACTION) {
/*
* load sensitivity index at this node point
* this routine determines the entry in the jacobian matrix which
* corresponds to this BC equation and this unknown - Jac_BC is a
* pointer to this entry *
* if it exists, add sens into matrix
*/
if (Dolphin[I][var] > 0) {
if (! doFullJac) {
ldof_var = ei->ln_to_first_dof[var][id];
if (ldof_var != -1) {
lec->J[ieqn][pvar][ldof_eqn][ldof_var] += penalty * d_func[var];
lec->J[ieqn][pvar][ldof_eqn][ldof_var] *= f_time;
}
} else {
if (doMeshMapping &&
(var == MESH_DISPLACEMENT1 || var == MESH_DISPLACEMENT2 || var == MESH_DISPLACEMENT3)) {
for (j = 0; j < ei->dof[var]; j++)
{
jk = dof_map[j];
lec->J[ieqn][pvar][ldof_eqn][jk] += penalty * d_kfunc[0][var][j];
lec->J[ieqn][pvar][ldof_eqn][jk] *= f_time;
}
} else {
for (j = 0; j < ei->dof[var]; j++)
{
lec->J[ieqn][pvar][ldof_eqn][j] += penalty * d_kfunc[0][var][j];
lec->J[ieqn][pvar][ldof_eqn][j] *= f_time;
}
}
}
} else {
/*
* if variable is not defined at this node, loop
* over all dof for this variable in this element
*/
for (j = 0; j < ei->dof[var]; j++) {
phi_j = bf[var]->phi[j];
lec->J[ieqn][pvar] [ldof_eqn][j] += penalty * d_func[var] * phi_j;
lec->J[ieqn][pvar] [ldof_eqn][j] *= f_time;
}
}
} else {
for (w = 0; w < pd->Num_Species_Eqn; w++) {
pvar = MAX_PROB_VAR + w;
if (Dolphin[I][var] > 0) {
ldof_var = ei->ln_to_first_dof[var][id];
if (ldof_var != -1) {
lec->J[ieqn][pvar] [ldof_eqn][ldof_var] += penalty * d_func[MAX_VARIABLE_TYPES + w];
lec->J[ieqn][pvar] [ldof_eqn][ldof_var] *= f_time;
}
}
/* if variable is not defined at this node,
* loop over all dof in this element */
else {
for (j = 0; j < ei->dof[var]; j++) {
phi_j = bf[var]->phi[j];
lec->J[ieqn][pvar] [ldof_eqn][j] += penalty * d_func[MAX_VARIABLE_TYPES + w] * phi_j;
lec->J[ieqn][pvar] [ldof_eqn][j] *= f_time;
}
}
} /* end of loop over species */
} /* end of if MASS_FRACTION */
} /* end of variable exists and condition is sensitive to it */
} /* end of loop over variable types */
} /* end of NEWTON */
} /* if (ldof_eqn != -1) */
} /* END of if (Res_BC != NULL), i.e. (index_eqn != -1) */
} /* END of if COLLOCATED BC */
/*****************************************************************************/
} /* END for (ibc = 0; (int) elem_side_bc->BC_input_id[ibc] != ...*/
/*****************************************************************************/
} /* END if (I < num_owned_nodes) */
/*****************************************************************************/
} /* END for (i = 0; i < (int) elem_side_bc->num_nodes_on_side; i++) */
/*****************************************************************************/
return(status);
} /* end of routine apply_point_colloc_bc() */
/*****************************************************************************/
/*****************************************************************************/
/*****************************************************************************/
void
moving_plane( int ielem_dim,
double *func,
double d_func[],
dbl *aa,
double time )
{
fplane ( ielem_dim, func, d_func, aa );
}
/*****************************************************************************/
/*****************************************************************************/
/*****************************************************************************/
void
fplane (int ielem_dim,
double *func,
double d_func[], /* dimensioned [MAX_VARIABLE_TYPES+MAX_CONC] */
dbl *aa) /* function parameters from data card */
{
/**************************** EXECUTION BEGINS *******************************/
int i;
if(af->Assemble_LSA_Mass_Matrix)
return;
*func = (double) aa[3];
for (i=0;i<ielem_dim;i++)
{
d_func[MESH_DISPLACEMENT1+i] = aa[i];
*func += (double) aa[i]*fv->x[i];
}
} /* END of routine fplane */
/*****************************************************************************/
#ifdef USE_CGM
void
sm_fplane (int ielem_dim,
double *func,
double d_func[], /* dimensioned [MAX_VARIABLE_TYPES+MAX_CONC] */
PlaneHandle *pHdl) /* Handle to a CGM Plane object */
{
EH(-1,"CGM bc_colloc.c:643 not implemented yet");
return;
} /* END of routine sm_fplane */
#endif
/*****************************************************************************/
void
fvelocity_profile (int var_flag,
int ielem_dim,
int velo_condition,
double *func,
double d_func[], /* defined [MAX_VARIABLE_TYPES + MAX_CONC] */
double p[], /* p[] are the parameters passed in through the input deck*/
double time) /* time at which bc's are evaluated */
{
if(af->Assemble_LSA_Mass_Matrix)
d_func[var_flag] = 0.0;
else
d_func[var_flag] = -1.0;
if( pd->e[R_MESH1] )
{
d_func[MESH_DISPLACEMENT1] =
dvelo_vary_fnc_d1(velo_condition, fv->x[0], fv->x[1], fv->x[2], p, time);
d_func[MESH_DISPLACEMENT2] =
dvelo_vary_fnc_d2(velo_condition, fv->x[0], fv->x[1], fv->x[2], p, time);
if (ielem_dim == 3) d_func[MESH_DISPLACEMENT3] =
dvelo_vary_fnc_d3(velo_condition, fv->x[0], fv->x[1], fv->x[2], p, time);
}
*func = velo_vary_fnc(velo_condition, fv->x[0], fv->x[1], fv->x[2], p, time);
*func -= fv->v[var_flag-VELOCITY1];
} /* END of routine fvelocity_profile */
/*****************************************************************************/
void
fspline (int ielem_dim,
double *func,
double d_func[], /* dimensioned [MAX_VARIABLE_TYPES + MAX_CONC] */
double p[], /* parameters to parameterize temperature eqn model*/
double time) /* time at which bc's are evaluated */
{
if(af->Assemble_LSA_Mass_Matrix)
return;
d_func[MESH_DISPLACEMENT1] =
dfncd1(fv->x[0], fv->x[1], fv->x[2], p, time);
d_func[MESH_DISPLACEMENT2] =
dfncd2(fv->x[0], fv->x[1], fv->x[2], p, time);
if (ielem_dim == 3) d_func[MESH_DISPLACEMENT3] =
dfncd3(fv->x[0], fv->x[1], fv->x[2], p, time);
*func = fnc(fv->x[0], fv->x[1], fv->x[2], p, time);
} /* END of routine fspline */
/*****************************************************************************/
void
fspline_rs (int ielem_dim,
double *func,
double d_func[], /* dimensioned [MAX_VARIABLE_TYPES + MAX_CONC] */
double p[], /* parameters to parameterize temperature eqn model*/
double time) /* time at which bc's are evaluated */
{
if(af->Assemble_LSA_Mass_Matrix)
return;
d_func[SOLID_DISPLACEMENT1] =
dfncd1(fv->x[0], fv->x[1], fv->x[2], p, time);
d_func[SOLID_DISPLACEMENT2] =
dfncd2(fv->x[0], fv->x[1], fv->x[2], p, time);
if (ielem_dim == 3) d_func[SOLID_DISPLACEMENT3] =
dfncd3(fv->x[0], fv->x[1], fv->x[2], p, time);
*func = fnc(fv->x[0], fv->x[1], fv->x[2], p, time);
} /* END of routine fspline_rs */
/*****************************************************************************/
void
fTmelting(double *func,
double d_func[], /* dimensioned MAX_VARIABLE_TYPES + MAX_CONC */
double a1) /* function parameter from data card */
{
if (af->Assemble_LSA_Mass_Matrix) return;
d_func[TEMPERATURE] = 1.;
*func = fv->T - a1;
} /* END of routine fTmelting */
/*****************************************************************************/
/******************************************************************************
fgeneralized_dirichlet() - simple pointwise collocation contribution for BCs
Function which adds simple contributions to boundary conditions.
This is used as a pointwise collocation boundary condition, and has
flexibility built in so that the condition can be applied to any equation
and be sensitive with respect to any unknown
Author: Rich Cairncross (1511)
Date: 22 December 1994
Revised:
******************************************************************************/
/*
* There are inconsistencies in the prototype declaration and how this function
* is invoked in bc_curve.c (circa line 1170).
*
* Based on other boundary conditions, I believe the declaration needs to be
* more like:
*
* double func[],
* double d_func[][MAX_VARIABLE_TYPES+MAX_CONC][MDE],
* ...
* Resolve this at some point. -PAS
*/
int
fgeneralized_dirichlet(double *func,
double d_func[], /* MAX_VARIABLE_TYPES + MAX_CONC */
const int gd_condition, /* denoting which condition
* applied */
const int bc_input_id,
const double tt, /* parameter to vary time integration
* from explicit (tt = 1) to
* implicit (tt = 0) */
const double dt) /* current time step size */
{
int jvar, wspec;
int index_var; /* Column index into the global stiffness matrix*/
dbl x_var; /* value of variable at this node */
dbl d_x_var; /* sensitivity of variable to nodal unknown */
dbl slope; /* slope of interpolated function in table */
dbl x_var_mp[1]; /* dummy variable for table lookup subroutines */
if(af->Assemble_LSA_Mass_Matrix)
return 0;
/* ---- Find variable number and species number */
jvar = BC_Types[bc_input_id].BC_Data_Int[2];
wspec = BC_Types[bc_input_id].BC_Data_Int[3];
/* put value of variable in GD Condition into x_var and put it's sensitivity in d_x_var */
index_var = load_variable( &x_var, &d_x_var, jvar, wspec, tt, dt);
/* Now add in contributions to residual vector and jacobian matrix */
switch(gd_condition)
{
case(GD_CONST_BC): /* x - c0 */
*func = (x_var - BC_Types[bc_input_id].BC_Data_Float[0] );
if (af->Assemble_Jacobian) {
d_func[index_var] = d_x_var;
}
break;
case(GD_LINEAR_BC): /* C1 x + c0 */
*func = (x_var* BC_Types[bc_input_id].BC_Data_Float[1]
+ BC_Types[bc_input_id].BC_Data_Float[0] );
if (af->Assemble_Jacobian) {
d_func[index_var] = d_x_var * BC_Types[bc_input_id].BC_Data_Float[1] ;
}
break;
case(GD_INVERSE_BC): /* C1/x + c0 */
*func = (BC_Types[bc_input_id].BC_Data_Float[1] / x_var
+ BC_Types[bc_input_id].BC_Data_Float[0] );
if (af->Assemble_Jacobian) {
d_func[index_var] = -d_x_var * BC_Types[bc_input_id].BC_Data_Float[1]/(x_var*x_var) ;
}
break;
case(GD_PARAB_BC): /* C2 x^2 + C1 x + c0 */
*func = (x_var * x_var * BC_Types[bc_input_id].BC_Data_Float[2]
+ x_var * BC_Types[bc_input_id].BC_Data_Float[1]
+ BC_Types[bc_input_id].BC_Data_Float[0] );
if (af->Assemble_Jacobian) {
d_func[index_var] = d_x_var * ( BC_Types[bc_input_id].BC_Data_Float[1]
+ 2. * x_var * BC_Types[bc_input_id].BC_Data_Float[2] );
}
break;
case(GD_PARAB_OFFSET_BC): /* C2 (x - C3)^2 + C1 (x - C3) + c0 */
*func = ((x_var - BC_Types[bc_input_id].BC_Data_Float[3])
* (x_var - BC_Types[bc_input_id].BC_Data_Float[3])
* BC_Types[bc_input_id].BC_Data_Float[2]
+ (x_var - BC_Types[bc_input_id].BC_Data_Float[3])
* BC_Types[bc_input_id].BC_Data_Float[1]
+ BC_Types[bc_input_id].BC_Data_Float[0] );
if (af->Assemble_Jacobian) {
d_func[index_var] = d_x_var * ( BC_Types[bc_input_id].BC_Data_Float[1]
+ 2. * (x_var - BC_Types[bc_input_id].BC_Data_Float[3])
* BC_Types[bc_input_id].BC_Data_Float[2] );
}
break;
case(GD_CIRC_BC): /* C2 ( x - C1 )^2 - c0^2 */
/* C2 represents ellipticity and C1 represents origin */
/* C0 is radius and should enter only one of the BC's */
*func = ( BC_Types[bc_input_id].BC_Data_Float[2] *
( x_var - BC_Types[bc_input_id].BC_Data_Float[1]) *
( x_var - BC_Types[bc_input_id].BC_Data_Float[1])
- BC_Types[bc_input_id].BC_Data_Float[0] * BC_Types[bc_input_id].BC_Data_Float[0]);
if (af->Assemble_Jacobian) {
d_func[index_var] = d_x_var * ( 2. * BC_Types[bc_input_id].BC_Data_Float[2] * (
x_var - BC_Types[bc_input_id].BC_Data_Float[1] ) );
}
break;
case(GD_POLYN_BC): /* up to 6th order polynomial */
/* C6 x^6 + C5 x^5 + C4 x^4 + C3 x^3 + C2 x^2 + C1 x + C0 */
*func = (x_var * x_var * x_var * x_var * x_var * x_var * BC_Types[bc_input_id].BC_Data_Float[6]
+ x_var * x_var * x_var * x_var * x_var * BC_Types[bc_input_id].BC_Data_Float[5]
+ x_var * x_var * x_var * x_var * BC_Types[bc_input_id].BC_Data_Float[4]
+ x_var * x_var * x_var * BC_Types[bc_input_id].BC_Data_Float[3]
+ x_var * x_var * BC_Types[bc_input_id].BC_Data_Float[2]
+ x_var * BC_Types[bc_input_id].BC_Data_Float[1]
+ BC_Types[bc_input_id].BC_Data_Float[0] );
/* printf("POLYN fit X,F = %f %f\n", x_var, *func); */
if (af->Assemble_Jacobian) {
d_func[index_var] = d_x_var * ( BC_Types[bc_input_id].BC_Data_Float[1]
+ 2. * x_var * BC_Types[bc_input_id].BC_Data_Float[2]
+ 3. * x_var * x_var * BC_Types[bc_input_id].BC_Data_Float[3]
+ 4. * x_var * x_var * x_var * BC_Types[bc_input_id].BC_Data_Float[4]
+ 5. * x_var * x_var * x_var * x_var * BC_Types[bc_input_id].BC_Data_Float[5]
+ 6. * x_var * x_var * x_var * x_var * x_var * BC_Types[bc_input_id].BC_Data_Float[6] );
}
break;
case(GD_TABLE_BC):
*func = BC_Types[bc_input_id].BC_Data_Float[0];
x_var_mp[0]=x_var;
*func *= interpolate_table( BC_Types[bc_input_id].table, x_var_mp, &slope, NULL );
if (af->Assemble_Jacobian)
{
d_func[index_var] = BC_Types[bc_input_id].BC_Data_Float[0]*slope*d_x_var;
}
break;
default:
return(-1);
}