forked from 2012ZGZYY/Dual_error_DG
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathdual_solver.cc
1347 lines (1195 loc) · 64.1 KB
/
dual_solver.cc
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 <deal.II/numerics/vector_tools.h>
#include <deal.II/numerics/solution_transfer.h>
#include <deal.II/numerics/matrix_tools.h>
#include <deal.II/lac/sparsity_pattern.h>
#include <deal.II/lac/trilinos_sparse_matrix.h>
#include <deal.II/lac/trilinos_vector.h>
#include <deal.II/lac/trilinos_precondition.h>
#include <deal.II/lac/trilinos_solver.h>
//#include <deal.II/lac/petsc_precondition.h>
//#include <deal.II/lac/petsc_solver.h>
#include <deal.II/lac/solver_gmres.h>
#include <deal.II/lac/sparse_direct.h>
#include <deal.II/lac/precondition_block.h>
#include <deal.II/lac/precondition.h> //I add this
#include <deal.II/fe/fe_tools.h>
#include "dual_solver.h"
#include <Sacado.hpp>
#include <iostream>
#include <fstream>
#include <vector>
#include <memory>
#include "equation.h"
using namespace dealii;
template<int dim>
DualSolver<dim>::
DualSolver(const Triangulation<dim>& triangulation,
const dealii::Mapping<dim>& mapping,
const FESystem<dim>& fe,
const Vector<double>& primal_solution,
const Parameters::AllParameters<dim>& parameters,
//const ConservationLaw<dim>& cons,
const DualFunctional::DualFunctionalBase<dim>& dual_functional,
const unsigned int& timestep_no
//const std::string& base_mesh
//const double& present_time,
//const double& end_time
)
:
mapping(mapping),
fe(fe),
parameters(parameters),
dof_handler(triangulation),
primal_solution_origin(primal_solution),
fe_degree_dual(fe.degree + 1), //TODO:
fe_dual(FE_DGQArbitraryNodes<dim>(QGauss<1>(fe_degree_dual+1)),
EulerEquations<dim>::n_components),
dof_handler_dual (triangulation),
quadrature_formula (fe_degree_dual+1),
face_quadrature_formula(fe_degree_dual+1),
//cons (cons),
triangulation(&triangulation), //triangulation is a pointer
dual_functional(&dual_functional),
timestep_no(timestep_no) //used in output_results()
//mpi_communicator(MPI_COMM_WORLD),
//this_mpi_process(Utilities::MPI::this_mpi_process(mpi_communicator))
//base_mesh(base_mesh),
//present_time(present_time),
//end_time(end_time)
{
//note, in the above constructor, the fe_dual must be initialized with FE_DGQArbitraryNodes, which
//is the same as in primal_solver.
}
template<int dim>
DualSolver<dim>::~DualSolver()
{
dof_handler_dual.clear();
}
template<int dim>
DualSolver<dim>::Integrator::
Integrator(const DoFHandler<dim>& dof_handler)
:
dof_info(dof_handler)
{}
//prepare the system_matrix, system_rhs and solution of the dual_problem. set them to correct size.
template<int dim>
void
DualSolver<dim>::setup_system()
{
dof_handler.distribute_dofs(fe);
dof_handler_dual.distribute_dofs(fe_dual);
std::cout<<"number of degrees of freedom in dual problem: "
<<dof_handler_dual.n_dofs()
<<std::endl;
solution_dual.reinit(dof_handler_dual.n_dofs());
dual_weights.reinit(dof_handler_dual.n_dofs());
system_rhs_dual.reinit(dof_handler_dual.n_dofs());
//------------------------------------------------------------------------------------------------
//set cell index. used to compute dt and mu_shock
/*
unsigned int index = 0;
for(typename Triangulation<dim>::active_cell_iterator cell = triangulation.begin_active();
cell != triangulation.end();
++cell, ++index)
cell->set_user_index(index);
*/
DynamicSparsityPattern c_sparsity(dof_handler_dual.n_dofs());
DoFTools::make_flux_sparsity_pattern(dof_handler_dual, c_sparsity);
sparsity_pattern_origin.copy_from(c_sparsity);
system_matrix_origin.reinit(sparsity_pattern_origin); //system_matrix_dual is its transpose
system_matrix_dual.reinit(sparsity_pattern_dual); //this sparsity_pattern will be changed after we
//compute system_matrix_dual
//generate an identity sparse matrix
sparsity_pattern_identity.reinit(dof_handler_dual.n_dofs(),dof_handler_dual.n_dofs(), 1);
sparsity_pattern_identity.compress();
system_matrix_identity.reinit(sparsity_pattern_identity); //throw error
for(SparseMatrix<double>::size_type i =0; i<system_matrix_identity.m(); i++)
system_matrix_identity.set(i, i, 1);
//then interpolate the primal_solution to the dual_space, we will use it as linearization point
//dof_handler and dof_handler_dual must be subscribed to the same triangulation
primal_solution.reinit(dof_handler_dual.n_dofs()); //denote primal_solution in dual_space, this is a member of dual_solver.
FETools::interpolate(dof_handler,
primal_solution_origin, //input
dof_handler_dual,
primal_solution); //output
}
template<int dim>
void
DualSolver<dim>::assemble_matrix()
{
//firstly, we initialize the integrator, this is almost the same as in ConservationLaw,
//the only difference is: we only need to assmble the system_matrix_dual here
//and don't need to worry about the system_rhs_dual.
Integrator integrator(dof_handler_dual);
//all the following quadrature formulas has n_gauss_points = fe_degree_dual+1;
integrator.info_box.initialize_gauss_quadrature(quadrature_formula.size(),
face_quadrature_formula.size(),
face_quadrature_formula.size());
integrator.info_box.initialize_update_flags ();
integrator.info_box.add_update_flags_all (update_values |
update_quadrature_points |
update_JxW_values);
integrator.info_box.add_update_flags_cell (update_gradients);
integrator.info_box.add_update_flags_boundary (update_normal_vectors | update_gradients);
integrator.info_box.add_update_flags_face (update_normal_vectors | update_gradients);
integrator.info_box.initialize (fe_dual, mapping);
//integrator.assembler.initialize (system_matrix_dual);
integrator.assembler.initialize (system_matrix_origin);
//call the worker functions on cells/boundaries/faces to assemble the dual_matrix
//system_matrix_dual = 0;
system_matrix_origin = 0;
std::cout<<"start to assemble the system_matrix_origin..."<<std::endl;
MeshWorker::loop<dim, dim, MeshWorker::DoFInfo<dim>, MeshWorker::IntegrationInfoBox<dim>>
(dof_handler_dual.begin_active(),
dof_handler_dual.end(),
integrator.dof_info,
integrator.info_box,
boost::bind(&DualSolver<dim>::integrate_cell_term, this, _1, _2), //TODO
boost::bind(&DualSolver<dim>::integrate_boundary_term, this, _1, _2),
boost::bind(&DualSolver<dim>::integrate_face_term, this, _1, _2, _3, _4),
integrator.assembler);
//system_matrix_dual.compress();
//=====================================================
std::cout<<"finish assembling the system_matrix_origin"<<std::endl;
// print out the system_matrix_origin(the transpose of the dual matrix)
/*
static int count = 0;
if(count == 0){
std::ofstream out("dual matrix",std::ios::out);
out.precision(15);
system_matrix_origin.print(out);
out.close();
count++;
}
*/
//=====================================================
}
template<int dim>
void
DualSolver<dim>::assemble_rhs()
{
//assemble rhs of the dual_system. input: primal_solution, output: rhs of the dual_problem
//rhs is a member of the dual class. note here for the PointValueEvaluation, primal_solution is not needed.
dual_functional->assemble_rhs(dof_handler_dual,
mapping,
primal_solution,
system_rhs_dual);
}
template<int dim>
std::pair<unsigned int, double>
DualSolver<dim>::solve()
{
/*get the transpose of system_matrix_origin, store the result in system_matrix_dual
system_matrix_origin is with sparsity_pattern_origin(with dof_handler_dual),
system_matrix_dual is with sparsity_pattern_dual,
system_matrix_identity is with sparsity_pattern_identity
*/
std::cout<<"start to transpose the system_matrix_origin..."<<std::endl;
system_matrix_origin.Tmmult(system_matrix_dual, system_matrix_identity);
std::cout<<"finish transposing the system_matrix_origin"<<std::endl;
/*
//print out the system_matrix_origin
static int count = 0;
if(count == 0){
std::ofstream out("dual matrix",std::ios::out);
out.precision(15);
system_matrix_origin.print(out);
out.close();
count++;
}
//system_matrix_dual.print_pattern(std::cout);
// static int count = 0;
if(count == 1){
std::ofstream out("dual matrix",std::ios::out);
out.precision(15);
system_matrix_dual.print(out);
out.close();
count++;
}
*/
/*
std::cout<<"system_matrix_dual dimension is "
<<system_matrix_dual.m()<<"x"<<system_matrix_dual.n()
<<std::endl;
*/
//now system_matrix_dual, system_rhs_dual are all done, solve it to get solution_dual
//---------------------------------------------------------------------------------------------------
//direct solver, the result is very good,
//but too slow and doesn't work for very large matrix.
/*
std::cout<<"start to solve the dual system..."<<std::endl;
SparseDirectUMFPACK solver;
solver.initialize(system_matrix_dual);
solver.vmult(solution_dual, system_rhs_dual);
std::cout<<"finish solving the dual system"<<std::endl;
return std::pair<unsigned int, double>(1,0);
*/
//---------------------------------------------------------------------------------------------------
//---------------------------------------------------------------------------------------------------
std::cout<<"start to solve the dual system..."<<std::endl;
//used to determine whether the iteration should be continued,
//30 is the maximum number of iteration steps before failure, default tolerance
//to determine success of the iteration is 1e-10
SolverControl solver_control(50);
//set the number of temparary vectors to 50, i.e. do a restart every 48 iterations,
//set right_preconditioning = true, set use_default_residual = true
SolverGMRES<>::AdditionalData gmres_data(50);
SolverGMRES<> solver(solver_control, gmres_data);
PreconditionBlockSSOR<SparseMatrix<double>, float> preconditioner;
preconditioner.initialize(system_matrix_dual, fe_dual.dofs_per_cell);
/* //Jacobi preconditioner
PreconditionBlockJacobi<SparseMatrix<double>, float> preconditioner;
preconditioner.initialize(system_matrix_dual,
PreconditionJacobi<SparseMatrix<double>>::AdditionalData(.6));
*/
try{
solver.solve(system_matrix_dual, solution_dual, system_rhs_dual, preconditioner);
}
catch(...){
std::cout<<" *** NO convergence in gmres when solving dual_problem ***\n";
}
std::cout<<"finish solving the dual system"<<std::endl;
return std::pair<unsigned int, double>(solver_control.last_step(),
solver_control.last_value());
}
//output solution_dual
template<int dim>
void
DualSolver<dim>::output_results()
{
int dual_iter = timestep_no/parameters.refine_iter_step;
std::string filename = "dual_solution-" + Utilities::to_string<unsigned int>(dual_iter, 4) + ".plt"; //TODO: file type
std::ofstream output (filename.c_str()); //convert std::string(C++ style) to const char*(C style)
//output the dual_solution
DataOut<dim> data_out;
data_out.attach_dof_handler (dof_handler_dual);
data_out.add_data_vector (solution_dual,
EulerEquations<dim>::component_names ()
/*DataOut<dim>::type_dof_data*/);
//output the dual_weights, this will be quite similar to the dual_solution, but with some difference
std::vector<std::string> dual_weights_names;
dual_weights_names.push_back ("XMomentum_dual_weights");
dual_weights_names.push_back ("YMomentum_dual_weights");
if(dim==3)
dual_weights_names.push_back ("ZMomentum_dual_weights");
dual_weights_names.push_back ("Density_dual_weights");
dual_weights_names.push_back ("Energy_dual_weights");
data_out.add_data_vector (dual_weights,
dual_weights_names);
data_out.build_patches (mapping, fe_degree_dual);
std::cout << "Writing file of dual_solution" << std::endl;
data_out.write_tecplot (output);
}
template<int dim>
void
DualSolver<dim>::compute_DWR_indicators(Vector<double>& refinement_indicators,
double& estimated_error_obj)
{
Assert (refinement_indicators.size()==triangulation->n_active_cells(),
ExcDimensionMismatch(refinement_indicators.size(),
triangulation->n_global_active_cells()));
//solve the dual_problem to get the dual_weights
setup_system();
assemble_matrix();
//======================================================
/*
static int count = 0;
if(count == 0){
std::ofstream out("dual matrix",std::ios::out);
out.precision(15);
system_matrix_origin.print(out);
out.close();
count++;
}
*/
//======================================================
assemble_rhs();
std::pair<unsigned int, double> convergence = solve();
std::cout<<"finish solving dual_problem ... "
<<std::endl
<<convergence.first<<" "<<convergence.second
<<std::endl;
//then in the following, compute the error indicator using dual_weights and residual
//------------------------compute the dual_weights---------------------------------
//dual_weights: z-z_h. z is in dual_space, z_h is z's interpolation into primal_space.
//consider this: why we must use z-z_h to denote the dual_weights, rather using z itself?
//follow function will compute (Id-I_h)z for a given dof-function z, where I_h is the
//interpolation from fe1 to fe2, i.e. from dual_space to primal_space
FETools::interpolation_difference (dof_handler_dual, //dof1
solution_dual, //z
fe, //fe2, i.e. primal_space
dual_weights); //z_differnce
dual_weights = solution_dual; //TODO
std::cout<<"finish computing the dual_weights..."<<std::endl;
//output the dual_solution
output_results();
//-----------------------prepare all the data_structure----------------------------
FEValues<dim> fe_values(fe_dual, quadrature_formula,
update_values| //TODO: quadrature_formula
update_gradients|
update_quadrature_points|
update_JxW_values);
const unsigned int n_q_points = quadrature_formula.size();
const unsigned int n_components = EulerEquations<dim>::n_components;
const unsigned int dofs_per_cell = fe_dual.dofs_per_cell;
std::vector<types::global_dof_index> dof_indices(dofs_per_cell);
std::vector<types::global_dof_index> dof_indices_neighbor(dofs_per_cell);
//data structures used to integrate on cells
std::vector<Vector<double>> dual_weights_cell_values(n_q_points,
Vector<double>(n_components));
std::vector<Vector<double>> cell_residual(n_q_points,
Vector<double>(n_components));
//std::vector<Vector<double>> W(n_q_points, Vector<double>(n_components));
//std::vector<std::vector<Vector<double>>> grad_W(n_q_points, std::vector<Vector<double>>(dim, Vector<double>(n_components)));
Table<2, double> W(n_q_points, n_components);
Table<3, double> grad_W(n_q_points, n_components, dim);
//note: why we use array rather than FullMatrix?
//sice FullMatrix need arguments to allocate memory each
//time such a matrix is created.
//std::vector<FullMatrix<double>> A_x(n_q_points, FullMatrix<double>(n_components,n_components)),
// A_y(n_q_points, FullMatrix<double>(n_components,n_components));
std::vector<std::array<std::array<double, n_components>, n_components>> A_x(n_q_points),
A_y(n_q_points);
std::vector<double> mu_shock_cell(n_q_points);
std::vector<double> residual_average(n_q_points);
std::vector<std::vector<Vector<double>>> grad_dual_weights(n_q_points, std::vector<Vector<double>>(dim, Vector<double>(n_components)));
/*std::vector<std::vector<Tensor<1, dim>>> grad_dual_weights(n_q_points);
for(unsigned int p=0; p!=n_q_points; ++p)
grad_dual_weights[p].resize(n_components); //denote d(z-z_h)/dx and d(z-z_h)/dy
*/
FEFaceValues<dim> fe_face_values_cell(fe_dual, face_quadrature_formula,
update_values|
update_quadrature_points|
update_gradients|
update_JxW_values|
update_normal_vectors),
fe_face_values_neighbor(fe_dual, face_quadrature_formula,
update_values|
update_quadrature_points|
update_gradients|
update_JxW_values|
update_normal_vectors);
FESubfaceValues<dim> fe_subface_values_cell(fe_dual, face_quadrature_formula,
update_values|
update_quadrature_points|
update_JxW_values|
update_normal_vectors),
fe_subface_values_neighbor(fe_dual,face_quadrature_formula,
update_values|
update_quadrature_points);
const unsigned int n_face_q_points = face_quadrature_formula.size();
//data structures used to integrate on faces
std::vector<Vector<double>> dual_weights_face_values(n_face_q_points, Vector<double>(n_components));
std::vector<Vector<double>> face_residual(n_face_q_points, Vector<double>(n_components));
//std::vector<Vector<double>> Wplus(n_face_q_points, Vector<double>(n_components));
//std::vector<Vector<double>> Wminus(n_face_q_points, Vector<double>(n_components));
Table<2, double> Wplus(n_face_q_points, n_components),
Wminus(n_face_q_points, n_components);
std::vector<Vector<double>> boundary_values(n_face_q_points,
Vector<double>(n_components));
//std::vector<Vector<double>> normal_flux(n_face_q_points, Vector<double>(n_components));
//std::vector<Vector<double>> numerical_flux(n_face_q_points, Vector<double>(n_components));
std::vector<std::array<double, n_components>> normal_flux(n_face_q_points); //should it be an array ?
std::vector<std::array<double, n_components>> numerical_flux(n_face_q_points);
//define a map from face_iterator to the face_term_error
//typename std::map<typename DoFHandler<dim>::face_iterator, Vector<double>> face_integrals;
typename DoFHandler<dim>::active_cell_iterator cell = dof_handler_dual.begin_active(),
endc = dof_handler_dual.end();
/*for(; cell!=endc; ++cell){
for(unsigned int face_no=0; face_no<GeometryInfo<dim>::faces_per_cell; ++face_no){
face_integrals[cell->face(face_no)].reinit(dim); //TODO: what mean?
face_integrals[cell->face(face_no)] = -1e20;
}
}*/
//----------------------------------loop-----------------------------------------
//cell = dof_handler_dual.begin_active();
unsigned int present_cell = 0;
for(; cell!=endc; ++cell, ++present_cell){
//-------------------------integrate_over_cell----------------------------
/*
std::cout<<"present_cell is "<<present_cell<<std::endl;
if(present_cell > 2694){
std::cout<<"the curent element is 2695"<<std::endl;
std::cout<<" --------------------------------- "<<std::endl;
}
*/
fe_values.reinit(cell);
cell->get_dof_indices (dof_indices); //store global_dof_indices corresponding to local_dof_indices
fe_values.get_function_values(dual_weights, dual_weights_cell_values);
//set initial values to 0 in every loop
for(unsigned int q_point=0; q_point<n_q_points; ++q_point){
residual_average[q_point] = 0;
for(unsigned int c=0; c<n_components; ++c){
W[q_point][c] = 0;
cell_residual[q_point][c] = 0;
for(unsigned int d=0; d<dim; d++)
grad_W[q_point][c][d] = 0; //dW/dx and dW/dy at quadrature_point, they are Vectors
}
for(unsigned int d=0; d<dim; d++)
grad_dual_weights[q_point][d] = 0;
}
//compute cell_residuals: df1(u)/dx+df2(u)/dy
for(unsigned int q_point=0; q_point<n_q_points; ++q_point){
//get grad_W and Jacobi matrix on all quadrature points
for(unsigned int i=0; i<dofs_per_cell; ++i){
const unsigned int c = fe_values.get_fe().system_to_component_index(i).first;
W[q_point][c] += primal_solution(dof_indices[i]) *
fe_values.shape_value_component(i,q_point,c);
for(unsigned int d=0; d<dim; d++)
grad_W[q_point][c][d] += primal_solution(dof_indices[i]) *
fe_values.shape_grad_component(i,q_point,c)[d];
}
EulerEquations<dim>::compute_Jacobi_matrix (W[q_point],
A_x[q_point],
A_y[q_point]); //get A_x, A_y
//use Jacobi matrix and grad_W to compute A_x*dW/dx+A_y*dW/dy, i.e. the residual
/*A_x[q_point].vmult_add(cell_residual[q_point], grad_W[q_point][0]); //cell_residual+=A_x*dW/dx+A_y*dW/dy
A_y[q_point].vmult_add(cell_residual[q_point], grad_W[q_point][1]);
*/
for(unsigned int row=0; row<n_components; row++){
for(unsigned int col=0; col<n_components; col++)
cell_residual[q_point][row] += A_x[q_point][row][col] * grad_W[q_point][col][0] +
A_y[q_point][row][col] * grad_W[q_point][col][1];
}
}
//cell_residuals are prepared, dual_weightes are prepared, now compute the integral
double cell_integral = 0;
for(unsigned int q_point=0; q_point<n_q_points; ++q_point){
cell_integral += cell_residual[q_point] *
dual_weights_cell_values[q_point] * //here is Vector*Vector
fe_values.JxW(q_point);
}
refinement_indicators(present_cell) -= cell_integral;
//compute mu_shock_cell, it's used to compute the integral term related to diffusion
for(unsigned int q_point=0; q_point<n_q_points; ++q_point){
for(unsigned int c=0; c<n_components; ++c)
residual_average[q_point] += std::abs(cell_residual[q_point][c]);
residual_average[q_point] /= n_components;
mu_shock_cell[q_point] = parameters.diffusion_coef*
std::pow(fe_values.get_cell()->diameter(), 2-parameters.diffusion_power)*
residual_average[q_point];
}
//integrate the cell term related to diffusion
for(unsigned int q_point=0; q_point<n_q_points; ++q_point){
for(unsigned int i=0; i<dofs_per_cell; ++i){
const unsigned int c = fe_values.get_fe().system_to_component_index(i).first;
for(unsigned int d=0; d<dim; d++)
grad_dual_weights[q_point][d][c] += dual_weights(dof_indices[i]) *
fe_values.shape_grad_component(i,q_point,c)[d];
}
}
//fe_values.get_function_gradients(dual_weights, grad_dual_weights);
cell_integral = 0;
for(unsigned int q_point=0; q_point<n_q_points; ++q_point){
for(unsigned int d=0; d<dim; ++d)
for(unsigned int c=0; c<n_components; ++c)
cell_integral += mu_shock_cell[q_point] *
grad_W[q_point][c][d] * grad_dual_weights[q_point][d][c] *
fe_values.JxW(q_point);
}
refinement_indicators(present_cell) -= cell_integral;
//-------------------------integrate_over_face----------------------------
for(unsigned int face_no=0; face_no<GeometryInfo<dim>::faces_per_cell; ++face_no){
//set data structure's initial value to 0
for(unsigned int q=0; q<n_face_q_points; ++q){
for(unsigned int c=0; c<n_components; ++c){
Wplus[q][c] = 0;
Wminus[q][c]= 0;
}
}
//if the face at boundary
if(cell->face(face_no)->at_boundary()){
const unsigned int& boundary_id = cell->face(face_no)->boundary_id();
//get Wplus at all quadrature points on current face
fe_face_values_cell.reinit(cell, face_no);
for(unsigned int q_point=0; q_point<n_face_q_points; ++q_point){
for(unsigned int i=0; i<dofs_per_cell; ++i){
const unsigned int c = fe_face_values_cell.get_fe().system_to_component_index(i).first;
Wplus[q_point][c] += primal_solution(dof_indices[i]) *
fe_face_values_cell.shape_value_component(i,q_point,c);
}
}
//fe_face_values_cell.get_function_values(primal_solution, Wplus);
//get Wminus at all quadrature points on current face
parameters.boundary_conditions[boundary_id].values.
vector_value_list(fe_face_values_cell.get_quadrature_points(),boundary_values);
typename EulerEquations<dim>::BoundaryKind boundary_kind =
parameters.boundary_conditions[boundary_id].kind;
for(unsigned int p=0; p<n_face_q_points; p++){
EulerEquations<dim>::compute_Wminus(boundary_kind,
fe_face_values_cell.normal_vector(p),
Wplus[p],
boundary_values[p],
Wminus[p]);
}
//now Wplus and Wminus are known, use them to compute normal_flux and numerical_normal_flux
for(unsigned int p=0; p<n_face_q_points; p++){
//compute the normal_flux and numerical_flux on current quadrature_point
EulerEquations<dim>::compute_normal_flux(Wplus[p],
fe_face_values_cell.normal_vector(p),
normal_flux[p]);
numerical_normal_flux(fe_face_values_cell.normal_vector(p),
Wplus[p],
Wminus[p],
numerical_flux[p]);
//compute the face_residual on current quadrature_point
for(unsigned int c=0; c<n_components; c++)
face_residual[p][c] = normal_flux[p][c] - numerical_flux[p][c];
}
//get the dual_weights
fe_face_values_cell.get_function_values(dual_weights, dual_weights_face_values);
//compute the face integral,i.e. the face_term_error
double face_integral = 0;
for(unsigned int q_point=0; q_point<n_face_q_points; ++q_point){
face_integral += face_residual[q_point] *
dual_weights_face_values[q_point] * //here is Vector*Vector
fe_face_values_cell.JxW(q_point);
}
refinement_indicators(present_cell) += face_integral;
}
//if the neighbor is finer(i.e. the current face has children)
else if(cell->face(face_no)->has_children()){
//
const unsigned int neighbor_neighbor = cell->neighbor_of_neighbor(face_no);
const typename DoFHandler<dim>::face_iterator face = cell->face(face_no);
for(unsigned int subface_no=0; subface_no<face->n_children(); ++subface_no){
//get Wplus. note that here Wplus only contains information on current subface
fe_subface_values_cell.reinit(cell, face_no, subface_no);
//fe_subface_values_cell.get_function_values(primal_solution, Wplus);
for(unsigned int q_point=0; q_point<n_face_q_points; ++q_point){
for(unsigned int i=0; i<dofs_per_cell; ++i){
const unsigned int c = fe_subface_values_cell.get_fe().
system_to_component_index(i).first;
Wplus[q_point][c] += primal_solution(dof_indices[i]) *
fe_subface_values_cell.shape_value_component(i,q_point,c);
}
}
//get Wminus
const typename DoFHandler<dim>::active_cell_iterator neighbor_child =
cell->neighbor_child_on_subface(face_no, subface_no);
Assert(neighbor_child->face(neighbor_neighbor)==cell->face(face_no)->child(subface_no),
ExcInternalError());
fe_face_values_neighbor.reinit(neighbor_child, neighbor_neighbor);
neighbor_child->get_dof_indices(dof_indices_neighbor);
//fe_face_values_neighbor.get_function_values(primal_solution, Wminus);
for(unsigned int q_point=0; q_point<n_face_q_points; ++q_point){
for(unsigned int i=0; i<dofs_per_cell; ++i){
const unsigned int c = fe_face_values_neighbor.get_fe().
system_to_component_index(i).first;
Wminus[q_point][c] += primal_solution(dof_indices_neighbor[i]) *
fe_face_values_neighbor.shape_value_component(i,q_point,c);
}
}
//use Wplus and Wminus on current subface to compute the flux and numerical flux
for(unsigned int p=0; p<n_face_q_points; p++){
//compute the normal_flux and numerical_flux on current quadrature_point
EulerEquations<dim>::compute_normal_flux(Wplus[p],
fe_subface_values_cell.normal_vector(p),
normal_flux[p]);
numerical_normal_flux(fe_subface_values_cell.normal_vector(p),
Wplus[p],
Wminus[p],
numerical_flux[p]);
//compute the face_residual on current quadrature_point
for(unsigned int c=0; c<n_components; c++)
face_residual[p][c] = normal_flux[p][c] - numerical_flux[p][c];
}
//get the dual_weights
fe_subface_values_cell.get_function_values(dual_weights, dual_weights_face_values);
//compute the face integral,i.e. the face_term_error
double subface_integral = 0;
for(unsigned int q_point=0; q_point<n_face_q_points; ++q_point){
subface_integral += face_residual[q_point] *
dual_weights_face_values[q_point] * //here is Vector*Vector
fe_subface_values_cell.JxW(q_point);
}
refinement_indicators(present_cell) += subface_integral;
}//end subface loop
}
//if the neighbor is as fine as our cell
else if(!cell->neighbor_is_coarser(face_no)){
//get Wplus
fe_face_values_cell.reinit(cell, face_no);
//fe_face_values_cell.get_function_values(primal_solution, Wplus);
for(unsigned int q_point=0; q_point<n_face_q_points; ++q_point){
for(unsigned int i=0; i<dofs_per_cell; ++i){
const unsigned int c = fe_face_values_cell.get_fe().
system_to_component_index(i).first;
Wplus[q_point][c] += primal_solution(dof_indices[i]) *
fe_face_values_cell.shape_value_component(i,q_point,c);
}
}
//get Wminus
Assert(cell->neighbor(face_no).state() == IteratorState::valid, ExcInternalError()); //check if the neighbor exists
const unsigned int neighbor_neighbor = cell->neighbor_of_neighbor(face_no);
const typename DoFHandler<dim>::active_cell_iterator neighbor = cell->neighbor(face_no);
fe_face_values_neighbor.reinit(neighbor, neighbor_neighbor);
neighbor->get_dof_indices(dof_indices_neighbor);
//fe_face_values_neighbor.get_function_values(primal_solution, Wminus);
for(unsigned int q_point=0; q_point<n_face_q_points; ++q_point){
for(unsigned int i=0; i<dofs_per_cell; ++i){
const unsigned int c = fe_face_values_neighbor.get_fe().
system_to_component_index(i).first;
Wminus[q_point][c] += primal_solution(dof_indices_neighbor[i]) *
fe_face_values_neighbor.shape_value_component(i,q_point,c);
}
}
//use Wplus and Wminus on current subface to compute the flux and numerical flux
for(unsigned int p=0; p<n_face_q_points; p++){
EulerEquations<dim>::compute_normal_flux(Wplus[p],
fe_face_values_cell.normal_vector(p),
normal_flux[p]);
numerical_normal_flux(fe_face_values_cell.normal_vector(p),
Wplus[p],
Wminus[p],
numerical_flux[p]);
//compute the face_residual on current quadrature_point
for(unsigned int c=0; c<n_components; c++)
face_residual[p][c] = normal_flux[p][c] - numerical_flux[p][c];
}
//get the dual_weights
fe_face_values_cell.get_function_values(dual_weights, dual_weights_face_values);
//compute the face integral,i.e. the face_term_error
double face_integral = 0;
for(unsigned int q_point=0; q_point<n_face_q_points; ++q_point){
face_integral += face_residual[q_point] *
dual_weights_face_values[q_point] * //here is Vector*Vector
fe_face_values_cell.JxW(q_point);
}
refinement_indicators(present_cell) += face_integral;
}
//if the neighbor is coarser
else{
//get Wplus
fe_face_values_cell.reinit(cell, face_no);
//fe_face_values_cell.get_function_values(primal_solution, Wplus);
for(unsigned int q_point=0; q_point<n_face_q_points; ++q_point){
for(unsigned int i=0; i<dofs_per_cell; ++i){
const unsigned int c = fe_face_values_cell.get_fe().
system_to_component_index(i).first;
Wplus[q_point][c] += primal_solution(dof_indices[i]) *
fe_face_values_cell.shape_value_component(i,q_point,c);
}
}
Assert(cell->neighbor(face_no).state() == IteratorState::valid, ExcInternalError());
Assert(cell->neighbor_is_coarser(face_no), ExcInternalError());
//get Wminus
const unsigned int neighbor_neighbor = cell->neighbor_of_coarser_neighbor(face_no).first;
const unsigned int subface_neighbor = cell->neighbor_of_coarser_neighbor(face_no).second;
const typename DoFHandler<dim>::active_cell_iterator neighbor = cell->neighbor(face_no);
fe_subface_values_neighbor.reinit(neighbor, neighbor_neighbor, subface_neighbor);
neighbor->get_dof_indices(dof_indices_neighbor);
//fe_subface_values_neighbor.get_function_values(primal_solution, Wminus);
for(unsigned int q_point=0; q_point<n_face_q_points; ++q_point){
for(unsigned int i=0; i<dofs_per_cell; ++i){
const unsigned int c = fe_subface_values_neighbor.get_fe().
system_to_component_index(i).first;
Wminus[q_point][c] += primal_solution(dof_indices_neighbor[i]) *
fe_subface_values_neighbor.shape_value_component(i,q_point,c);
}
}
//use Wplus and Wminus on current subface to compute the flux and numerical flux
for(unsigned int p=0; p<n_face_q_points; p++){
EulerEquations<dim>::compute_normal_flux(Wplus[p],
fe_face_values_cell.normal_vector(p),
normal_flux[p]);
numerical_normal_flux(fe_face_values_cell.normal_vector(p),
Wplus[p],
Wminus[p],
numerical_flux[p]);
//compute the face_residual on current quadrature_point
for(unsigned int c=0; c<n_components; c++)
face_residual[p][c] = normal_flux[p][c] - numerical_flux[p][c]; //here is Vector-Vector
}
//get the dual_weights
fe_face_values_cell.get_function_values(dual_weights, dual_weights_face_values);
//compute the face integral,i.e. the face_term_error
double face_integral = 0;
for(unsigned int q_point=0; q_point<n_face_q_points; ++q_point){
face_integral += face_residual[q_point] *
dual_weights_face_values[q_point] * //here is Vector*Vector
fe_face_values_cell.JxW(q_point);
}
refinement_indicators(present_cell) += face_integral;
}
}//end face_loop
}// end cell_loop
std::cout<<"finish computing the error indicator"<<std::endl;
//following is the estimated error in the objective functional.
estimated_error_obj = std::accumulate(refinement_indicators.begin(),
refinement_indicators.end(), 0.);
std::cout<<"estimated error = "
<<estimated_error_obj
<<std::endl;
//set the indicators on each cell to positive
for(Vector<double>::iterator i = refinement_indicators.begin(); i!=refinement_indicators.end(); ++i)
*i = std::fabs(*i);
}
//
template<int dim>
void
DualSolver<dim>::integrate_cell_term(DoFInfo& dinfo, CellInfo& info)
{
FullMatrix<double>& local_matrix = dinfo.matrix(0).matrix;
//Vector<double>& local_vector = dinfo.vector(0).block(0);
std::vector<unsigned int>& dof_indices = dinfo.indices;
const FEValuesBase<dim>& fe_v = info.fe_values();
const unsigned int dofs_per_cell = fe_v.dofs_per_cell;
const unsigned int n_q_points = fe_v.n_quadrature_points;
//const double time_step = dt(cell_number(dinfo.cell));
//Table<2,Sacado::Fad::DFad<double> >
// dWdt (n_q_points, EulerEquations<dim>::n_components);
Table<2, Sacado::Fad::DFad<double>>
W (n_q_points, EulerEquations<dim>::n_components);
Table<3, Sacado::Fad::DFad<double>>
grad_W (n_q_points, EulerEquations<dim>::n_components, dim);
// extract primal_solution and store as local variable. here primal_solution is used.
std::vector<Sacado::Fad::DFad<double> > independent_local_dof_values(dofs_per_cell);
for (unsigned int i=0; i<dofs_per_cell; ++i)
independent_local_dof_values[i] = primal_solution(dof_indices[i]);
for (unsigned int i=0; i<dofs_per_cell; ++i)
independent_local_dof_values[i].diff (i, dofs_per_cell);
for (unsigned int q=0; q<n_q_points; ++q)
for (unsigned int i=0; i<dofs_per_cell; ++i)
{
const unsigned int c = fe_v.get_fe().system_to_component_index(i).first;
Sacado::Fad::DFad<double> dof = independent_local_dof_values[i];
W[q][c] += dof * fe_v.shape_value_component(i, q, c);
for (unsigned int d = 0; d < dim; d++)
grad_W[q][c][d] += dof * fe_v.shape_grad_component(i, q, c)[d];
}
std::vector<std::array<std::array<Sacado::Fad::DFad<double>, dim>,
EulerEquations<dim>::n_components>>
flux(n_q_points);
std::vector<std::array<Sacado::Fad::DFad<double>, EulerEquations<dim>::n_components>>
forcing(n_q_points);
//--------------------------------------------------------------------------------------
//I add following variables in order to compute the diffusion terms
//--------------------------------------------------------------------------------------
//R denotes the residual vectors on quadrature points of current cell
Table<2,Sacado::Fad::DFad<double>> R(n_q_points, EulerEquations<dim>::n_components);
//get the jacobi on quadrature points
std::vector<std::array<std::array<Sacado::Fad::DFad<double>,
EulerEquations<dim>::n_components>,
EulerEquations<dim>::n_components>>
A_x(n_q_points),
A_y(n_q_points);
//--------------------------------------------------------------------------------------
for (unsigned int q=0; q<n_q_points; ++q)
{
EulerEquations<dim>::compute_flux_matrix (W[q], flux[q]);
EulerEquations<dim>::compute_forcing_vector (W[q], forcing[q]);
EulerEquations<dim>::compute_Jacobi_matrix (W[q], A_x[q], A_y[q]); //get A_x & A_y
}
// Get current cell number
//const unsigned int cell_no = cell_number (dinfo.cell);
//--------------------------------------------------------------------------------------
// I add this vector to denote the diffusion coefficient
std::vector<Sacado::Fad::DFad<double>> mu_shock_cell(n_q_points);
std::vector<Sacado::Fad::DFad<double>> R_average(n_q_points);
for(unsigned int point=0; point<n_q_points; ++point){
for(unsigned int c=0; c<EulerEquations<dim>::n_components; c++){
//R[point][c] += dWdt[point][c];
for(unsigned int j=0; j<EulerEquations<dim>::n_components; j++){ //the j-th component of the variable vector W
R[point][c] += A_x[point][c][j]*grad_W[point][j][0];
R[point][c] += A_y[point][c][j]*grad_W[point][j][1];
}
}
}//finish computing Residual vector R(point, components)
for(unsigned int point=0; point<n_q_points; ++point){
for(unsigned int c=0; c<EulerEquations<dim>::n_components; ++c){
R_average[point] += R[point][c]*R[point][c]; //L2 norm
}
R_average[point] = std::sqrt(R_average[point]);
mu_shock_cell[point] = parameters.diffusion_coef*
std::pow(fe_v.get_cell()->diameter(),
2-parameters.diffusion_power)*
R_average[point];
}
/*
//compute the average of mu_shock_cell on a cell just for visualization
Sacado::Fad::DFad<double> mu_shock_ave = 0;
for(unsigned int point=0; point<n_q_points; ++point)
mu_shock_ave += mu_shock_cell[point];
mu_shock_ave /= n_q_points;
mu_shock(cell_no) = mu_shock_ave.val(); */
//--------------------------------------------------------------------------------------
for (unsigned int i=0; i<dofs_per_cell; ++i)
{
Sacado::Fad::DFad<double> F_i = 0;
const unsigned int
component_i = fe_v.get_fe().system_to_component_index(i).first;
// The residual for each row (i) will be accumulating
// into this fad variable. At the end of the assembly
// for this row, we will query for the sensitivities
// to this variable and add them into the Jacobian.
for (unsigned int point=0; point<n_q_points; ++point)
{
for (unsigned int d=0; d<dim; d++)
F_i -= flux[point][component_i][d] *
fe_v.shape_grad_component(i, point, component_i)[d] *
fe_v.JxW(point);
// Diffusion term for shocks
/*
if(parameters.diffusion_coef > 0.0)
for (unsigned int d=0; d<dim; d++)
F_i += mu_shock(cell_no) *
grad_W[point][component_i][d] *
fe_v.shape_grad_component(i, point, component_i)[d] *
fe_v.JxW(point);
*/
//my new way of adding diffusiton terms
if(parameters.diffusion_coef > 0.0)
for (unsigned int d=0; d<dim; d++)
F_i += mu_shock_cell[point] *
grad_W[point][component_i][d] *
fe_v.shape_grad_component(i, point, component_i)[d] *
fe_v.JxW(point);
/*
F_i -= parameters.gravity *
forcing[point][component_i] *
fe_v.shape_value_component(i, point, component_i) *
fe_v.JxW(point);
*/
}
for (unsigned int k=0; k<dofs_per_cell; ++k)
local_matrix (i, k) += F_i.fastAccessDx(k);
}
}
template<int dim>
void
DualSolver<dim>::integrate_boundary_term(DoFInfo& dinfo, CellInfo& info)
{
FullMatrix<double>& local_matrix = dinfo.matrix(0).matrix;
//Vector<double>& local_vector = dinfo.vector(0).block(0);
std::vector<unsigned int>& dof_indices = dinfo.indices;
const unsigned int& face_no = dinfo.face_number;
const double& face_diameter = dinfo.face->diameter();
const unsigned int& boundary_id = dinfo.face->boundary_id();
const FEValuesBase<dim>& fe_v = info.fe_values();
const unsigned int n_q_points = fe_v.n_quadrature_points;
const unsigned int dofs_per_cell = fe_v.dofs_per_cell;
std::vector<Sacado::Fad::DFad<double> >
independent_local_dof_values (dofs_per_cell);
const unsigned int n_independent_variables = dofs_per_cell;
// Get cell number of two cells adjacent to this face
//const unsigned int cell_no = cell_number (dinfo.cell);
for (unsigned int i = 0; i < dofs_per_cell; i++)
{
independent_local_dof_values[i] = primal_solution(dof_indices[i]);
independent_local_dof_values[i].diff(i, n_independent_variables);
}
// Conservative variable value at face
Table<2,Sacado::Fad::DFad<double> >
Wplus (n_q_points, EulerEquations<dim>::n_components),
Wminus (n_q_points, EulerEquations<dim>::n_components);
/* TODO:ADIFF
Table<3,Sacado::Fad::DFad<double> >
grad_Wplus (n_q_points, EulerEquations<dim>::n_components, dim);
*/
for (unsigned int q=0; q<n_q_points; ++q)
for (unsigned int i=0; i<dofs_per_cell; ++i)
{
const unsigned int c = fe_v.get_fe().system_to_component_index(i).first;
Sacado::Fad::DFad<double> dof = independent_local_dof_values[i];
Wplus[q][c] += dof * fe_v.shape_value_component(i, q, c);
/* TODO:ADIFF
for (unsigned int d = 0; d < dim; d++)
grad_Wplus[q][c][d] += dof *
fe_v.shape_grad_component(i, q, c)[d];
*/