-
Notifications
You must be signed in to change notification settings - Fork 0
/
residue.cpp
1733 lines (1374 loc) · 59.9 KB
/
residue.cpp
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
#include "common.hpp"
template <typename Derived>
void getProjectorMatrix(MatrixBase<Derived> & P, int n_nodes, int const* nodes, Vec const& Vec_x_, double t, AppCtx const& app)
{
int const dim = app.dim;
Mesh const* mesh = &*app.mesh;
//DofHandler const* dof_handler = &*app.dof_handler;
std::vector<int> const& dirichlet_tags = app.dirichlet_tags;
//std::vector<int> const& neumann_tags = app.neumann_tags ;
//std::vector<int> const& interface_tags = app.interface_tags;
std::vector<int> const& solid_tags = app.solid_tags ;
std::vector<int> const& triple_tags = app.triple_tags ;
//std::vector<int> const& periodic_tags = app.periodic_tags ;
std::vector<int> const& feature_tags = app.feature_tags ;
//Vec const& Vec_normal = app.Vec_normal;
P.setIdentity();
Tensor I(dim,dim);
Tensor Z(dim,dim);
Vector X(dim);
Vector normal(dim);
int dofs[dim];
int tag;
Point const* point;
I.setIdentity();
Z.setZero();
// NODES
for (int i = 0; i < n_nodes; ++i)
{
point = mesh->getNodePtr(nodes[i]);
tag = point->getTag();
//m = point->getPosition() - mesh->numVerticesPerCell();
//cell = mesh->getCellPtr(point->getIncidCell());
if (is_in(tag,feature_tags))
{
app.getNodeDofs(&*point, DH_MESH, VAR_M, dofs);
VecGetValues(Vec_x_, dim, dofs, X.data());
P.block(i*dim,i*dim,dim,dim) = feature_proj(X,t,tag);
continue;
}
else
if (is_in(tag,solid_tags) || is_in(tag, triple_tags))
{
app.getNodeDofs(&*point, DH_MESH, VAR_M, dofs);
VecGetValues(Vec_x_, dim, dofs, X.data());
normal = solid_normal(X,t,tag);
P.block(i*dim,i*dim,dim,dim) = I - normal*normal.transpose();
//P.block(i*dim,i*dim,dim,dim) = Z;
}
else
if (is_in(tag,dirichlet_tags))
{
P.block(i*dim,i*dim,dim,dim) = Z;
}
} // end nodes
}
// ******************************************************************************
// FORM FUNCTION
// ******************************************************************************
PetscErrorCode AppCtx::formFunction(SNES /*snes*/, Vec Vec_up_k, Vec Vec_fun)
{
//PetscErrorCode ierr;
int null_space_press_dof=-1;
int iter;
SNESGetIterationNumber(snes,&iter);
if (!iter)
{
converged_times=0;
}
if (force_pressure && (iter<2))
{
Vector X(dim);
Vector X_new(dim);
if (behaviors & BH_Press_grad_elim)
{
cell_iterator cell = mesh->cellBegin();
dof_handler[DH_UNKS].getVariable(VAR_P).getCellDofs(&null_space_press_dof, &*cell);
// fix the initial guess
VecSetValue(Vec_up_k, null_space_press_dof, 0.0, INSERT_VALUES);
}
else
{
point_iterator point = mesh->pointBegin();
while (!mesh->isVertex(&*point))
++point;
int x_dofs[3];
dof_handler[DH_MESH].getVariable(VAR_M).getVertexDofs(x_dofs, &*point);
VecGetValues(Vec_x_1, dim, x_dofs, X_new.data());
VecGetValues(Vec_x_0, dim, x_dofs, X.data());
X = .5*(X+X_new);
dof_handler[DH_UNKS].getVariable(VAR_P).getVertexDofs(&null_space_press_dof, &*point);
// fix the initial guess
VecSetValue(Vec_up_k, null_space_press_dof, pressure_exact(X,current_time+.5*dt,point->getTag()), INSERT_VALUES);
}
Assembly(Vec_up_k);
}
// checking:
if (null_space_press_dof < 0 && force_pressure==1 && (iter<2))
{
cout << "force_pressure: somthing is wrong ..." << endl;
throw;
}
Mat *JJ = &Mat_Jac;
//PetscErrorCode ierr;
VecZeroEntries(Vec_fun);
MatZeroEntries(*JJ);
// LOOP NAS CÉLULAS
#ifdef FEP_HAS_OPENMP
FEP_PRAGMA_OMP(parallel default(none) shared(Vec_up_k,Vec_fun,cout,null_space_press_dof,JJ))
#endif
{
VectorXd FUloc(n_dofs_u_per_cell);
VectorXd FPloc(n_dofs_p_per_cell);
/* local data */
int tag;
MatrixXd u_coefs_c_mid_trans(dim, n_dofs_u_per_cell/dim); // n+utheta
MatrixXd u_coefs_c_old(n_dofs_u_per_cell/dim, dim); // n
MatrixXd u_coefs_c_old_trans(dim,n_dofs_u_per_cell/dim); // n
MatrixXd u_coefs_c_new(n_dofs_u_per_cell/dim, dim); // n+1
MatrixXd u_coefs_c_new_trans(dim,n_dofs_u_per_cell/dim); // n+1
MatrixXd v_coefs_c_mid(nodes_per_cell, dim); // mesh velocity; n
MatrixXd v_coefs_c_mid_trans(dim,nodes_per_cell); // mesh velocity; n
VectorXd p_coefs_c_new(n_dofs_p_per_cell); // n+1
VectorXd p_coefs_c_old(n_dofs_p_per_cell); // n
VectorXd p_coefs_c_mid(n_dofs_p_per_cell); // n
//MatrixXd x_coefs_c_mid(nodes_per_cell, dim); // n+utheta
MatrixXd x_coefs_c_mid_trans(dim, nodes_per_cell); // n+utheta
MatrixXd x_coefs_c_new(nodes_per_cell, dim); // n+1
MatrixXd x_coefs_c_new_trans(dim, nodes_per_cell); // n+1
MatrixXd x_coefs_c_old(nodes_per_cell, dim); // n
MatrixXd x_coefs_c_old_trans(dim, nodes_per_cell); // n
Tensor F_c_mid(dim,dim); // n+utheta
Tensor invF_c_mid(dim,dim); // n+utheta
Tensor invFT_c_mid(dim,dim); // n+utheta
Tensor F_c_old(dim,dim); // n
Tensor invF_c_old(dim,dim); // n
Tensor invFT_c_old(dim,dim); // n
Tensor F_c_new(dim,dim); // n+1
Tensor invF_c_new(dim,dim); // n+1
Tensor invFT_c_new(dim,dim); // n+1
//Tensor F_c_new(dim,dim); // n+1
//Tensor invF_c_new(dim,dim); // n+1
//Tensor invFT_c_new(dim,dim); // n+1
//Tensor F_c_old(dim,dim); // n
//Tensor invF_c_old(dim,dim); // n
//Tensor invFT_c_old(dim,dim); // n
/* All variables are in (n+utheta) by default */
MatrixXd dxphi_c(n_dofs_u_per_cell/dim, dim);
MatrixXd dxphi_c_new(dxphi_c);
MatrixXd dxpsi_c(n_dofs_p_per_cell, dim);
MatrixXd dxpsi_c_new(dxpsi_c);
MatrixXd dxqsi_c(nodes_per_cell, dim);
Vector dxbble(dim);
Vector dxbble_new(dim);
Tensor dxU(dim,dim); // grad u
Tensor dxU_old(dim,dim); // grad u
Tensor dxU_new(dim,dim); // grad u
Tensor dxUb(dim,dim); // grad u bble
Vector dxP_new(dim); // grad p
Vector Xqp(dim);
Vector Xc(dim); // cell center; to compute CR element
Vector Uqp(dim);
Vector Ubqp(dim); // bble
Vector Uqp_old(dim); // n
Vector Uqp_new(dim); // n+1
Vector Vqp(dim);
Vector Uconv_qp(dim);
Vector dUdt(dim);
double Pqp_new;
double Pqp;
double bble_integ=0;
//VectorXd FUloc(n_dofs_u_per_cell); // subvetor da função f (parte de U)
//VectorXd FPloc(n_dofs_p_per_cell); // subvetor da função f (parte de P)
VectorXi cell_nodes(nodes_per_cell);
double J_mid;
double J_new, J_old;
double JxW_mid;
double JxW_new, JxW_old;
double weight;
double visc=-1; // viscosity
double cell_volume;
double hk2;
double tauk=0;
double delk=0;
double delta_cd;
double rho;
MatrixXd Aloc(n_dofs_u_per_cell, n_dofs_u_per_cell);
MatrixXd Gloc(n_dofs_u_per_cell, n_dofs_p_per_cell);
MatrixXd Dloc(n_dofs_p_per_cell, n_dofs_u_per_cell);
MatrixXd Eloc(n_dofs_p_per_cell, n_dofs_p_per_cell); // GSL, BC
MatrixXd Cloc(n_dofs_u_per_cell, n_dofs_p_per_cell); // GSL
Tensor iBbb(dim, dim); // BC, i : inverse ..it is not the inverse to CR element
MatrixXd Bbn(dim, n_dofs_u_per_cell); // BC
MatrixXd Bnb(n_dofs_u_per_cell, dim); // BC
MatrixXd Dpb(n_dofs_p_per_cell, dim); // BC
MatrixXd Gbp(dim, n_dofs_p_per_cell); // BC
MatrixXd Gnx(n_dofs_u_per_cell, dim); // CR ;; suffix x means p gradient
Vector FUb(dim); // BC
Vector FPx(dim); // pressure gradient
Vector force_at_mid(dim);
Vector Res(dim); // residue
Tensor dResdu(dim,dim); // residue derivative
Tensor const I(Tensor::Identity(dim,dim));
Vector vec(dim); // temp
Tensor Ten(dim,dim); // temp
VectorXi mapU_c(n_dofs_u_per_cell);
VectorXi mapU_r(n_dofs_u_per_corner);
VectorXi mapP_c(n_dofs_p_per_cell);
VectorXi mapP_r(n_dofs_p_per_corner);
// mesh velocity
VectorXi mapM_c(dim*nodes_per_cell);
VectorXi mapM_f(dim*nodes_per_facet);
VectorXi mapM_r(dim*nodes_per_corner);
MatrixXd Prj(n_dofs_u_per_cell,n_dofs_u_per_cell); // projector matrix
//VectorXi cell_nodes(nodes_per_cell);
const int tid = omp_get_thread_num();
const int nthreads = omp_get_num_threads();
cell_iterator cell = mesh->cellBegin(tid,nthreads);
cell_iterator cell_end = mesh->cellEnd(tid,nthreads);
//cell_iterator cell = mesh->cellBegin();
//cell_iterator cell_end = mesh->cellEnd();
for (; cell != cell_end; ++cell)
{
tag = cell->getTag();
// mapeamento do local para o global:
//
dof_handler[DH_MESH].getVariable(VAR_M).getCellDofs(mapM_c.data(), &*cell);
dof_handler[DH_UNKS].getVariable(VAR_U).getCellDofs(mapU_c.data(), &*cell);
dof_handler[DH_UNKS].getVariable(VAR_P).getCellDofs(mapP_c.data(), &*cell);
VecGetValues(Vec_v_mid, mapM_c.size(), mapM_c.data(), v_coefs_c_mid.data());
VecGetValues(Vec_x_0, mapM_c.size(), mapM_c.data(), x_coefs_c_old.data());
VecGetValues(Vec_x_1, mapM_c.size(), mapM_c.data(), x_coefs_c_new.data());
VecGetValues(Vec_up_0, mapU_c.size(), mapU_c.data(), u_coefs_c_old.data());
VecGetValues(Vec_up_k, mapU_c.size(), mapU_c.data(), u_coefs_c_new.data());
VecGetValues(Vec_up_k, mapP_c.size(), mapP_c.data(), p_coefs_c_new.data());
VecGetValues(Vec_up_0, mapP_c.size(), mapP_c.data(), p_coefs_c_old.data());
// get nodal coordinates of the old and new cell
mesh->getCellNodesId(&*cell, cell_nodes.data());
//mesh->getNodesCoords(cell_nodes.begin(), cell_nodes.end(), x_coefs_c.data());
//x_coefs_c_trans = x_coefs_c_mid_trans;
v_coefs_c_mid_trans = v_coefs_c_mid.transpose();
x_coefs_c_old_trans = x_coefs_c_old.transpose();
x_coefs_c_new_trans = x_coefs_c_new.transpose();
u_coefs_c_old_trans = u_coefs_c_old.transpose();
u_coefs_c_new_trans = u_coefs_c_new.transpose();
u_coefs_c_mid_trans = utheta*u_coefs_c_new_trans + (1.-utheta)*u_coefs_c_old_trans;
x_coefs_c_mid_trans = utheta*x_coefs_c_new_trans + (1.-utheta)*x_coefs_c_old_trans;
p_coefs_c_mid = utheta*p_coefs_c_new + (1.-utheta)*p_coefs_c_old;
visc = muu(tag);
rho = pho(Xqp,tag);
Aloc.setZero();
Gloc.setZero();
Dloc.setZero();
FUloc.setZero();
FPloc.setZero();
Eloc.setZero();
if (behaviors & BH_bble_condens_PnPn) // reset matrices
{
iBbb.setZero();
Bnb.setZero();
Gbp.setZero();
FUb.setZero();
Bbn.setZero();
Dpb.setZero();
}
if(behaviors & BH_GLS)
{
cell_volume = 0;
for (int qp = 0; qp < n_qpts_cell; ++qp) {
F_c_mid = x_coefs_c_mid_trans * dLqsi_c[qp];
J_mid = determinant(F_c_mid,dim);
cell_volume += J_mid * quadr_cell->weight(qp);
}
hk2 = cell_volume / pi; // element size
double const uconv = (u_coefs_c_old - v_coefs_c_mid).lpNorm<Infinity>();
tauk = 4.*visc/hk2 + 2.*rho*uconv/sqrt(hk2);
tauk = 1./tauk;
if (dim==3)
tauk *= 0.1;
delk = 4.*visc + 2.*rho*uconv*sqrt(hk2);
//delk = 0;
Eloc.setZero();
Cloc.setZero();
}
if (behaviors & BH_bble_condens_CR)
{
bble_integ = 0;
Gnx.setZero();
iBbb.setZero();
Bnb.setZero();
FUb.setZero();
FPx.setZero();
Bbn.setZero();
cell_volume = 0;
Xc.setZero();
for (int qp = 0; qp < n_qpts_cell; ++qp) {
F_c_mid = x_coefs_c_mid_trans * dLqsi_c[qp];
J_mid = determinant(F_c_mid,dim);
Xqp = x_coefs_c_mid_trans * qsi_c[qp];
cell_volume += J_mid * quadr_cell->weight(qp);
Xc += J_mid * quadr_cell->weight(qp) * Xqp;
}
Xc /= cell_volume;
}
// Quadrature
for (int qp = 0; qp < n_qpts_cell; ++qp)
{
//F_c_new = x_coefs_c_new_trans * dLqsi_c[qp];
//inverseAndDet(F_c_new,dim,invF_c_new,J_new);
//invFT_c_new= invF_c_new.transpose();
//F_c_old = x_coefs_c_old_trans * dLqsi_c[qp];
//inverseAndDet(F_c_old,dim,invF_c_old,J_old);
//invFT_c_old= invF_c_old.transpose();
F_c_mid = x_coefs_c_mid_trans * dLqsi_c[qp];
F_c_old = x_coefs_c_old_trans * dLqsi_c[qp];
F_c_new = x_coefs_c_new_trans * dLqsi_c[qp];
inverseAndDet(F_c_mid,dim,invF_c_mid,J_mid);
inverseAndDet(F_c_old,dim,invF_c_old,J_old);
inverseAndDet(F_c_new,dim,invF_c_new,J_new);
invFT_c_mid= invF_c_mid.transpose();
invFT_c_old= invF_c_old.transpose();
invFT_c_new= invF_c_new.transpose();
dxphi_c_new = dLphi_c[qp] * invF_c_new;
dxphi_c = dLphi_c[qp] * invF_c_mid;
dxpsi_c_new = dLpsi_c[qp] * invF_c_new;
dxpsi_c = dLpsi_c[qp] * invF_c_mid;
dxqsi_c = dLqsi_c[qp] * invF_c_mid;
dxP_new = dxpsi_c.transpose() * p_coefs_c_new;
dxU = u_coefs_c_mid_trans * dxphi_c; // n+utheta
dxU_old = u_coefs_c_old_trans * dLphi_c[qp] * invF_c_old; // n
dxU_new = u_coefs_c_new_trans * dLphi_c[qp] * invF_c_new; // n+1
Xqp = x_coefs_c_mid_trans * qsi_c[qp]; // coordenada espacial (x,y,z) do ponto de quadratura
Uqp = u_coefs_c_mid_trans * phi_c[qp]; //n+utheta
Uqp_new = u_coefs_c_new_trans * phi_c[qp]; //n+utheta
Uqp_old = u_coefs_c_old_trans * phi_c[qp]; //n+utheta
Pqp_new = p_coefs_c_new.dot(psi_c[qp]);
Pqp = p_coefs_c_mid.dot(psi_c[qp]);
Vqp = v_coefs_c_mid_trans * qsi_c[qp];
Uconv_qp = Uqp - Vqp;
//Uconv_qp = Uqp_old;
dUdt = (Uqp_new-Uqp_old)/dt;
force_at_mid = force(Xqp,current_time+utheta*dt,tag);
weight = quadr_cell->weight(qp);
JxW_mid = J_mid*weight;
JxW_old = J_old*weight;
JxW_new = J_new*weight;
//~ if (mesh->getCellId(&*cell) == 0)
//~ {
//~ printf("cHEcKKKKKKKKKKK!!\n");
//~ cout << "x coefs mid:" << endl;
//~ cout << x_coefs_c_mid_trans.transpose() << endl;
//~ }
if (J_mid < 1.e-14)
{
FEP_PRAGMA_OMP(critical)
//if (tid==0)
{
printf("in formCellFunction:\n");
std::cout << "erro: jacobiana da integral não invertível: ";
std::cout << "J_mid = " << J_mid << endl;
cout << "trans matrix:\n" << F_c_mid << endl;
cout << "x coefs mid:" << endl;
cout << x_coefs_c_mid_trans.transpose() << endl;
cout << "-----" << endl;
cout << "cell id: " << mesh->getCellId(&*cell) << endl;
cout << "cell Contig id: " << mesh->getCellContigId( mesh->getCellId(&*cell) ) << endl;
cout << "cell nodes:\n" << cell_nodes.transpose() << endl;
cout << "mapM :\n" << mapM_c.transpose() << endl;
throw;
}
}
for (int i = 0; i < n_dofs_u_per_cell/dim; ++i)
{
for (int c = 0; c < dim; ++c)
{
FUloc(i*dim + c) += JxW_mid*
( rho*(unsteady*dUdt(c) + has_convec*Uconv_qp.dot(dxU.row(c)))*phi_c[qp][i] + // aceleração
visc*dxphi_c.row(i).dot(dxU.row(c) + dxU.col(c).transpose()) //rigidez
) -
JxW_mid*force_at_mid(c)*phi_c[qp][i] -// força +
//JxW_new*Pqp_new*dxphi_c_new(i,c); // pressão
//JxW_mid*Pqp*dxphi_c(i,c); // pressão
JxW_mid*Pqp_new*dxphi_c(i,c); // pressão
for (int j = 0; j < n_dofs_u_per_cell/dim; ++j)
{
for (int d = 0; d < dim; ++d)
{
delta_cd = c==d;
Aloc(i*dim + c, j*dim + d) += JxW_mid*
( has_convec*phi_c[qp][i]*utheta *rho*( delta_cd*Uconv_qp.dot(dxphi_c.row(j)) + dxU(c,d)*phi_c[qp][j] ) // advecção
+ unsteady* delta_cd*rho*phi_c[qp][i]*phi_c[qp][j]/dt // time derivative
+ utheta*visc*( delta_cd * dxphi_c.row(i).dot(dxphi_c.row(j)) + dxphi_c(i,d)*dxphi_c(j,c)) ); // rigidez
}
}
for (int j = 0; j < n_dofs_p_per_cell; ++j)
{
Gloc(i*dim + c,j) -= JxW_mid * psi_c[qp][j]* dxphi_c(i,c);
//Gloc(i*dim + c,j) -= utheta*JxW_mid * psi_c[qp][j]* dxphi_c(i,c);
//Gloc(i*dim + c,j) -= JxW_new * psi_c[qp][j]* dxphi_c_new(i,c);
//Dloc(j, i*dim + c) -= JxW_new * psi_c[qp][j]* dxphi_c_new(i,c);
Dloc(j, i*dim + c) -= utheta*JxW_mid * psi_c[qp][j]* dxphi_c(i,c);
}
}
//FUloc.segment(i*dim,dim) += JxW_mid*
// ( phi_c[qp][i]*rho* ( dUdt + has_convec* dxU*Uconv_qp ) // aceleração
// + visc* ( dxU + dxU.transpose() )*dxphi_c.row(i).transpose() //rigidez
// - Pqp_new*dxphi_c.row(i).transpose() // pressão
// - phi_c[qp][i]* force_at_mid ); // força
}
for (int i = 0; i < n_dofs_p_per_cell; ++i)
FPloc(i) -= JxW_mid* dxU.trace()*psi_c[qp][i];
//FPloc(i) -= JxW_new* dxU_new.trace() *psi_c[qp][i];
// ----------------
//
// STABILIZATION
//
// ----------------
if (behaviors & (BH_bble_condens_PnPn | BH_bble_condens_CR))
{
dxbble = invFT_c_mid * dLbble[qp];
for (int c = 0; c < dim; c++)
{
for (int j = 0; j < n_dofs_u_per_cell/dim; j++)
{
for (int d = 0; d < dim; d++)
{
delta_cd = c==d;
Bbn(c, j*dim + d) += JxW_mid*
( has_convec*bble[qp]*utheta *rho*( delta_cd*Uconv_qp.dot(dxphi_c.row(j)) + dxU(c,d)*phi_c[qp][j] ) // convective
+ unsteady*delta_cd*rho*bble[qp]*phi_c[qp][j]/dt // time derivative
+ utheta*visc*(delta_cd * dxphi_c.row(j).dot(dxbble) + dxphi_c(j,c)*dxbble(d)) ); // rigidez
Bnb(j*dim + d, c) += JxW_mid*
( has_convec*phi_c[qp][j]*utheta *rho*( delta_cd*Uconv_qp.dot(dxbble) ) // convective
+ delta_cd*rho*phi_c[qp][j]*bble[qp]/dt * unsteady // time derivative
+ utheta*visc*(delta_cd * dxphi_c.row(j).dot(dxbble) + dxphi_c(j,c)*dxbble(d)) ); // rigidez
}
}
if (behaviors & BH_bble_condens_PnPn)
for (int j = 0; j < n_dofs_p_per_cell; ++j)
Gbp(c, j) -= JxW_mid*psi_c[qp][j]*dxbble(c);
}
for (int c = 0; c < dim; c++)
{
for (int d = 0; d < dim; d++)
{
delta_cd = c==d;
iBbb(c, d) += JxW_mid*
( has_convec*bble[qp]*utheta *rho*( delta_cd*Uconv_qp.dot(dxbble) ) // convective
+ delta_cd*rho*bble[qp]*bble[qp]/dt * unsteady // time derivative
+ utheta*visc*(delta_cd* dxbble.dot(dxbble) + dxbble(d)*dxbble(c)) ); // rigidez
}
FUb(c) += JxW_mid*
( bble[qp]*rho*(dUdt(c)*unsteady + has_convec*Uconv_qp.dot(dxU.row(c))) + // time derivative + convective
visc*dxbble.dot(dxU.row(c) + dxU.col(c).transpose()) - //rigidez
Pqp_new*dxbble(c) - // pressão
force_at_mid(c)*bble[qp] ); // força
}
}
else
if(behaviors & BH_GLS)
{
Res = rho*( dUdt * unsteady + has_convec*dxU*Uconv_qp) + dxP_new - force_at_mid;
for (int j = 0; j < n_dofs_u_per_cell/dim; ++j)
{
dResdu = unsteady*(rho*phi_c[qp][j]/dt)*I + has_convec*rho*utheta*( phi_c[qp][j]*dxU + Uconv_qp.dot(dxphi_c.row(j))*I );
for (int i = 0; i < n_dofs_p_per_cell; ++i)
{
vec = dxpsi_c.row(i).transpose();
vec = dResdu.transpose()*vec;
vec = -JxW_mid*tauk* vec;
for (int d = 0; d < dim; d++)
Dloc(i, j*dim + d) += vec(d);
// atençao nos indices
vec = JxW_mid*tauk* has_convec*rho*Uconv_qp.dot(dxphi_c.row(j))* dxpsi_c.row(i).transpose();
for (int d = 0; d < dim; d++)
Cloc(j*dim + d,i) += vec(d);
}
for (int i = 0; i < n_dofs_u_per_cell/dim; ++i)
{
// supg term
Ten = JxW_mid*tauk* has_convec*( utheta*rho*phi_c[qp][j]*Res*dxphi_c.row(i) + rho*Uconv_qp.dot(dxphi_c.row(i))*dResdu );
// divergence term
Ten+= JxW_mid*delk*utheta*dxphi_c.row(i).transpose()*dxphi_c.row(j);
for (int c = 0; c < dim; ++c)
for (int d = 0; d < dim; ++d)
Aloc(i*dim + c, j*dim + d) += Ten(c,d);
}
}
for (int i = 0; i < n_dofs_p_per_cell; ++i)
for (int j = 0; j < n_dofs_p_per_cell; ++j)
Eloc(i,j) -= tauk*JxW_mid * dxphi_c.row(i).dot(dxphi_c.row(j));
for (int i = 0; i < n_dofs_u_per_cell/dim; i++)
{
vec = JxW_mid*( has_convec*tauk* rho* Uconv_qp.dot(dxphi_c.row(i)) * Res + delk*dxU.trace()*dxphi_c.row(i).transpose() );
for (int c = 0; c < dim; c++)
FUloc(i*dim + c) += vec(c);
}
for (int i = 0; i < n_dofs_p_per_cell; ++i)
FPloc(i) -= JxW_mid *tauk* dxpsi_c.row(i).dot(Res);
//FPloc(i) -= JxW_mid *tauk* dxpsi_c.row(i).dot(dxP_new - force_at_mid); // somente laplaciano da pressao
}
if (behaviors & BH_bble_condens_CR)
{
bble_integ += JxW_mid*bble[qp];
for (int c = 0; c < dim; ++c)
{
for (int i = 0; i < n_dofs_u_per_cell/dim; ++i)
for (int j = 0; j < dim; ++j) // pressure gradient
Gnx(i*dim + c,j) -= JxW_mid* (Xqp(j) - Xc(j))*dxphi_c(i,c);
FPx(c) -= JxW_mid* dxU.trace()*(Xqp(c) - Xc(c));
}
}
} // fim quadratura
//Dloc += utheta*Gloc.transpose();
//
// stabilization
//
if (behaviors & BH_bble_condens_PnPn)
{
//iBbb = iBbb.inverse().eval();
invert(iBbb,dim);
FUloc = FUloc - Bnb*iBbb*FUb;
FPloc = FPloc - utheta*Gbp.transpose()*iBbb*FUb;
Dpb = utheta*Gbp.transpose();
// correções com os coeficientes da bolha
Ubqp = -utheta*iBbb*FUb; // U bolha no tempo n+utheta
for (int qp = 0; qp < n_qpts_cell; ++qp)
{
F_c_mid = x_coefs_c_mid_trans * dLqsi_c[qp];
inverseAndDet(F_c_mid, dim, invF_c_mid,J_mid);
invFT_c_mid= invF_c_mid.transpose();
Uqp = u_coefs_c_mid_trans * phi_c[qp]; //n+utheta
dxbble = invFT_c_mid * dLbble[qp];
dxUb = Ubqp*dxbble.transpose();
weight = quadr_cell->weight(qp);
JxW_mid = J_mid*weight;
for (int j = 0; j < n_dofs_u_per_cell/dim; ++j)
{
for (int i = 0; i < n_dofs_u_per_cell/dim; ++i)
{
Ten = has_convec*JxW_mid*rho*utheta*phi_c[qp][i]* phi_c[qp][j] * dxUb; // advecção
for (int c = 0; c < dim; ++c)
for (int d = 0; d < dim; ++d)
Aloc(i*dim + c, j*dim + d) += Ten(c,d);
}
Ten = has_convec*JxW_mid* rho*utheta* bble[qp] * phi_c[qp][j] *dxUb; // advecção
for (int c = 0; c < dim; ++c)
for (int d = 0; d < dim; ++d)
Bbn(c, j*dim + d) += Ten(c,d);
}
} // fim quadratura 2 vez
Aloc -= Bnb*iBbb*Bbn;
Gloc -= Bnb*iBbb*Gbp;
Dloc -= Dpb*iBbb*Bbn;
Eloc = -Dpb*iBbb*Gbp;
}
if(behaviors & BH_GLS)
{
Gloc += Cloc;
}
if(behaviors & BH_bble_condens_CR)
{
Ubqp.setZero();
//for (int i = 0; i < Gnx.cols(); ++i)
// for (int j = 0; j < Gnx.rows(); ++j)
// Ubqp(i) += Gnx(j,i)*u_coefs_c_new(j);
//Ubqp /= -bble_integ;
//Ubqp *= utheta;
for (int c = 0; c < dim; ++c)
for (int i = 0; i < n_dofs_u_per_cell/dim; ++i)
for (int j = 0; j < dim; ++j) // pressure gradient
//Gnx(i*dim + c,j) -= JxW_mid* (Xqp(j) - Xc(j))*dxphi_c(i,c);
Ubqp(j) += Gnx(i*dim + c,j) * u_coefs_c_mid_trans(c,i);
Ubqp /= -bble_integ;
Ubqp *= utheta;
//Ubqp = -Gnx.transpose()*u_coefs_c_new; // U bolha no tempo n+utheta
for (int qp = 0; qp < n_qpts_cell; ++qp)
{
F_c_mid = x_coefs_c_mid_trans * dLqsi_c[qp];
inverseAndDet(F_c_mid, dim, invF_c_mid,J_mid);
invFT_c_mid= invF_c_mid.transpose();
Uqp = u_coefs_c_mid_trans * phi_c[qp]; //n+utheta
dxbble = invFT_c_mid * dLbble[qp];
dxUb = Ubqp*dxbble.transpose();
weight = quadr_cell->weight(qp);
JxW_mid = J_mid*weight;
for (int j = 0; j < n_dofs_u_per_cell/dim; ++j)
{
for (int i = 0; i < n_dofs_u_per_cell/dim; ++i)
{
Ten = has_convec*JxW_mid*rho*utheta*phi_c[qp][i]* phi_c[qp][j] * dxUb; // advecção
for (int c = 0; c < dim; ++c)
for (int d = 0; d < dim; ++d)
Aloc(i*dim + c, j*dim + d) += Ten(c,d);
}
Ten = has_convec*JxW_mid* rho*utheta* bble[qp] * phi_c[qp][j] *dxUb; // advecção
for (int c = 0; c < dim; ++c)
for (int d = 0; d < dim; ++d)
Bbn(c, j*dim + d) += Ten(c,d);
}
} // fim quadratura 2 vez
double const a = 1./(bble_integ*bble_integ);
double const b = 1./bble_integ;
Aloc += utheta*a*Gnx*iBbb*Gnx.transpose() - utheta*b*Bnb*Gnx.transpose() - b*Gnx*Bbn;
FUloc += a*Gnx*iBbb*FPx - b*Bnb*FPx - b*Gnx*FUb;
}
// Projection - to force non-penetrarion bc
mesh->getCellNodesId(&*cell, cell_nodes.data());
getProjectorMatrix(Prj, nodes_per_cell, cell_nodes.data(), Vec_x_1, current_time+dt, *this);
FUloc = Prj*FUloc;
Aloc = Prj*Aloc*Prj;
Gloc = Prj*Gloc;
Dloc = Dloc*Prj;
if (force_pressure)
{
for (int i = 0; i < mapP_c.size(); ++i)
{
if (mapP_c(i) == null_space_press_dof)
{
Gloc.col(i).setZero();
Dloc.row(i).setZero();
FPloc(i) = 0;
Eloc.col(i).setZero();
Eloc.row(i).setZero();
break;
}
}
}
#ifdef FEP_HAS_OPENMP
FEP_PRAGMA_OMP(critical)
#endif
{
VecSetValues(Vec_fun, mapU_c.size(), mapU_c.data(), FUloc.data(), ADD_VALUES);
VecSetValues(Vec_fun, mapP_c.size(), mapP_c.data(), FPloc.data(), ADD_VALUES);
MatSetValues(*JJ, mapU_c.size(), mapU_c.data(), mapU_c.size(), mapU_c.data(), Aloc.data(), ADD_VALUES);
MatSetValues(*JJ, mapU_c.size(), mapU_c.data(), mapP_c.size(), mapP_c.data(), Gloc.data(), ADD_VALUES);
MatSetValues(*JJ, mapP_c.size(), mapP_c.data(), mapU_c.size(), mapU_c.data(), Dloc.data(), ADD_VALUES);
MatSetValues(*JJ, mapP_c.size(), mapP_c.data(), mapP_c.size(), mapP_c.data(), Eloc.data(), ADD_VALUES);
}
}
}
// LOOP NAS FACES DO CONTORNO
//~ FEP_PRAGMA_OMP(parallel default(none) shared(Vec_up_k,Vec_fun,cout))
{
int tag;
bool is_neumann;
bool is_surface;
bool is_solid;
//MatrixXd u_coefs_c_new(n_dofs_u_per_facet/dim, dim);
//VectorXd p_coefs_f(n_dofs_p_per_facet);
MatrixXd u_coefs_f_mid_trans(dim, n_dofs_u_per_facet/dim); // n+utheta
MatrixXd u_coefs_f_old(n_dofs_u_per_facet/dim, dim); // n
MatrixXd u_coefs_f_new(n_dofs_u_per_facet/dim, dim); // n+1
MatrixXd u_coefs_f_old_trans(dim,n_dofs_u_per_facet/dim); // n
MatrixXd u_coefs_f_new_trans(dim,n_dofs_u_per_facet/dim); // n+1
MatrixXd x_coefs_f_mid_trans(dim, n_dofs_v_per_facet/dim); // n+utheta
MatrixXd x_coefs_f_new(n_dofs_v_per_facet/dim, dim); // n+1
MatrixXd x_coefs_f_new_trans(dim, n_dofs_v_per_facet/dim); // n+1
MatrixXd x_coefs_f_old(n_dofs_v_per_facet/dim, dim); // n
MatrixXd x_coefs_f_old_trans(dim, n_dofs_v_per_facet/dim); // n
MatrixXd noi_coefs_f_new(n_dofs_v_per_facet/dim, dim); // normal interpolada em n+1
MatrixXd noi_coefs_f_new_trans(dim, n_dofs_v_per_facet/dim); // normal interpolada em n+1
Tensor F_f_mid(dim,dim-1); // n+utheta
Tensor invF_f_mid(dim-1,dim); // n+utheta
Tensor fff_f_mid(dim-1,dim-1); // n+utheta; fff = first fundamental form
//Tensor invFT_c_mid(dim,dim); // n+utheta
MatrixXd Aloc_f(n_dofs_u_per_facet, n_dofs_u_per_facet);
VectorXd FUloc(n_dofs_u_per_facet);
MatrixXd tmp(n_dofs_u_per_facet,n_dofs_u_per_facet);
VectorXi mapU_f(n_dofs_u_per_facet);
VectorXi mapP_f(n_dofs_p_per_facet);
VectorXi mapM_f(dim*nodes_per_facet);
MatrixXd Prj(n_dofs_u_per_facet,n_dofs_u_per_facet);
MatrixXd dxphi_f(n_dofs_u_per_facet/dim, dim);
Tensor dxU_f(dim,dim); // grad u
Vector Xqp(dim);
Vector Xqp2(dim);
Vector Xqp_new(dim);
Vector Xqp_old(dim);
Vector Uqp(dim);
Vector Uqp_new(dim);
Vector Uqp_old(dim);
//VectorXd FUloc(n_dofs_u_per_facet);
VectorXi facet_nodes(nodes_per_facet);
Vector normal(dim);
Vector noi(dim); // normal interpolada
Vector some_vec(dim);
double J_mid=0,JxW_mid;
double weight=0;
//double visc;
//double rho;
Vector Uqp_solid(dim);
Vector traction_(dim);
//~ const int tid = omp_get_thread_num();
//~ const int nthreads = omp_get_num_threads();
//~
//~ facet_iterator facet = mesh->facetBegin(tid,nthreads);
//~ facet_iterator facet_end = mesh->facetEnd(tid,nthreads);
// LOOP NAS FACES DO CONTORNO
facet_iterator facet = mesh->facetBegin();
facet_iterator facet_end = mesh->facetEnd();
if (neumann_tags.size() != 0 || interface_tags.size() != 0 || solid_tags.size() != 0)
for (; facet != facet_end; ++facet)
{
tag = facet->getTag();
is_neumann = is_in(tag, neumann_tags);
is_surface = is_in(tag, interface_tags);
is_solid = is_in(tag, solid_tags);
if ((!is_neumann) && (!is_surface) && (!is_solid))
//PetscFunctionReturn(0);
continue;
// mapeamento do local para o global:
//
dof_handler[DH_UNKS].getVariable(VAR_U).getFacetDofs(mapU_f.data(), &*facet);
dof_handler[DH_UNKS].getVariable(VAR_P).getFacetDofs(mapP_f.data(), &*facet);
dof_handler[DH_MESH].getVariable(VAR_M).getFacetDofs(mapM_f.data(), &*facet);
VecGetValues(Vec_normal, mapM_f.size(), mapM_f.data(), noi_coefs_f_new.data());
VecGetValues(Vec_x_0, mapM_f.size(), mapM_f.data(), x_coefs_f_old.data());
VecGetValues(Vec_x_1, mapM_f.size(), mapM_f.data(), x_coefs_f_new.data());
VecGetValues(Vec_up_0, mapU_f.size(), mapU_f.data(), u_coefs_f_old.data());
VecGetValues(Vec_up_k , mapU_f.size(), mapU_f.data(), u_coefs_f_new.data());
// get nodal coordinates of the old and new cell
mesh->getFacetNodesId(&*facet, facet_nodes.data());
x_coefs_f_old_trans = x_coefs_f_old.transpose();
x_coefs_f_new_trans = x_coefs_f_new.transpose();
u_coefs_f_old_trans = u_coefs_f_old.transpose();
u_coefs_f_new_trans = u_coefs_f_new.transpose();
noi_coefs_f_new_trans = noi_coefs_f_new.transpose();
u_coefs_f_mid_trans = utheta*u_coefs_f_new_trans + (1.-utheta)*u_coefs_f_old_trans;
x_coefs_f_mid_trans = utheta*x_coefs_f_new_trans + (1.-utheta)*x_coefs_f_old_trans;
FUloc.setZero();
Aloc_f.setZero();
//visc = muu(tag);
//rho = pho(Xqp,tag);
//noi_coefs_f_new_trans = x_coefs_f_mid_trans;
for (int qp = 0; qp < n_qpts_facet; ++qp)
{
F_f_mid = x_coefs_f_mid_trans * dLqsi_f[qp];
if (dim==2)
{
normal(0) = +F_f_mid(1,0);
normal(1) = -F_f_mid(0,0);
normal.normalize();
}
else
{
normal = cross(F_f_mid.col(0), F_f_mid.col(1));
normal.normalize();
}
fff_f_mid.resize(dim-1,dim-1);
fff_f_mid = F_f_mid.transpose()*F_f_mid;
J_mid = sqrt(fff_f_mid.determinant());
invF_f_mid = fff_f_mid.inverse()*F_f_mid.transpose();
weight = quadr_facet->weight(qp);
JxW_mid = J_mid*weight;
Xqp = x_coefs_f_mid_trans * qsi_f[qp]; // coordenada espacial (x,y,z) do ponto de quadratura
dxphi_f = dLphi_f[qp] * invF_f_mid;
dxU_f = u_coefs_f_mid_trans * dxphi_f; // n+utheta
Uqp = u_coefs_f_mid_trans * phi_f[qp];
noi = noi_coefs_f_new_trans * qsi_f[qp];
if (is_neumann)
{
//Vector no(Xqp);
//no.normalize();
//traction_ = utheta*(traction(Xqp,current_time+dt,tag)) + (1.-utheta)*traction(Xqp,current_time,tag);
traction_ = traction(Xqp, normal, current_time + dt*utheta,tag);
//traction_ = (traction(Xqp,current_time,tag) +4.*traction(Xqp,current_time+dt/2.,tag) + traction(Xqp,current_time+dt,tag))/6.;
for (int i = 0; i < n_dofs_u_per_facet/dim; ++i)
{
for (int c = 0; c < dim; ++c)
{
FUloc(i*dim + c) -= JxW_mid * traction_(c) * phi_f[qp][i] ; // força
}
}
}
if (is_surface)
{
//Vector no(Xqp);
//no.normalize();
for (int i = 0; i < n_dofs_u_per_facet/dim; ++i)
{
for (int c = 0; c < dim; ++c)
{
//FUloc(i*dim + c) += JxW_mid *gama(Xqp,current_time,tag)*(dxphi_f(i,c) + 0*(unsteady*dt) *dxU_f.row(c).dot(dxphi_f.row(i))); // correto
FUloc(i*dim + c) += JxW_mid *gama(Xqp,current_time,tag)*dxphi_f(i,c);
//for (int d = 0; d < dim; ++d)
// FUloc(i*dim + c) += JxW_mid * gama(Xqp,current_time,tag)* ( (c==d?1:0) - noi(c)*noi(d) )* dxphi_f(i,d) ;
//FUloc(i*dim + c) += JxW_mid * gama(Xqp,current_time,tag)* ( unsteady*dt *dxU_f.row(c).dot(dxphi_f.row(i)));
}
}
if (false) // semi-implicit term
{
for (int i = 0; i < n_dofs_u_per_facet/dim; ++i)
for (int j = 0; j < n_dofs_u_per_facet/dim; ++j)
for (int c = 0; c < dim; ++c)
Aloc_f(i*dim + c, j*dim + c) += utheta*JxW_mid* (unsteady*dt) *gama(Xqp,current_time,tag)*dxphi_f.row(i).dot(dxphi_f.row(j));
}
}
if (is_solid)
{
Uqp_solid = solid_veloc(Xqp, current_time+utheta*dt, tag);
for (int i = 0; i < n_dofs_u_per_facet/dim; ++i)
{
for (int c = 0; c < dim; ++c)
{
FUloc(i*dim + c) += JxW_mid *beta_diss()*(Uqp(c)-Uqp_solid(c))*phi_f[qp][i];
//FUloc(i*dim + c) += x_coefs_f_old_trans.norm()*beta_diss()*Uqp(c)*phi_f[qp][i];
for (int j = 0; j < n_dofs_u_per_facet/dim; ++j)
Aloc_f(i*dim + c, j*dim + c) += utheta*JxW_mid *beta_diss()*phi_f[qp][j]*phi_f[qp][i];
}
}
}
} // end quadratura
// Projection - to force non-penetrarion bc
mesh->getFacetNodesId(&*facet, facet_nodes.data());
getProjectorMatrix(Prj, nodes_per_facet, facet_nodes.data(), Vec_x_1, current_time+dt, *this);
FUloc = Prj*FUloc;
Aloc_f = Prj*Aloc_f*Prj;
//~ FEP_PRAGMA_OMP(critical)
{
VecSetValues(Vec_fun, mapU_f.size(), mapU_f.data(), FUloc.data(), ADD_VALUES);
MatSetValues(*JJ, mapU_f.size(), mapU_f.data(), mapU_f.size(), mapU_f.data(), Aloc_f.data(), ADD_VALUES);
}
}
} // end parallel