MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
complexweakform.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#include "complexweakform.hpp"
13
14namespace mfem
15{
16
18{
19 trial_integs_r.SetSize(trial_fes.Size(), test_fecols.Size());
20 trial_integs_i.SetSize(trial_fes.Size(), test_fecols.Size());
21 for (int i = 0; i < trial_integs_r.NumRows(); i++)
22 {
23 for (int j = 0; j < trial_integs_r.NumCols(); j++)
24 {
27 }
28 }
29
30 test_integs_r.SetSize(test_fecols.Size(), test_fecols.Size());
31 test_integs_i.SetSize(test_fecols.Size(), test_fecols.Size());
32 for (int i = 0; i < test_integs_r.NumRows(); i++)
33 {
34 for (int j = 0; j < test_integs_r.NumCols(); j++)
35 {
38 }
39 }
40
41 lfis_r.SetSize(test_fecols.Size());
42 lfis_i.SetSize(test_fecols.Size());
43 for (int j = 0; j < lfis_r.Size(); j++)
44 {
47 }
48
50
51 mat_r = mat_e_r = NULL;
52 mat_i = mat_e_i = NULL;
55 width = height;
56
57 initialized = true;
58 static_cond = nullptr;
59
61 {
62 Bmat.SetSize(mesh->GetNE());
63 fvec.SetSize(mesh->GetNE());
64 }
65}
66
68{
71 dof_offsets[0] = 0;
72 tdof_offsets[0] = 0;
73 for (int i =0; i<nblocks; i++)
74 {
75 dof_offsets[i+1] = trial_fes[i]->GetVSize();
76 tdof_offsets[i+1] = trial_fes[i]->GetTrueVSize();
77 }
80}
81
82// Allocate SparseMatrix and RHS
84{
85 if (static_cond) { return; }
86
88 mat_r->owns_blocks = 1;
90 mat_i->owns_blocks = 1;
91
92 for (int i = 0; i < mat_r->NumRowBlocks(); i++)
93 {
94 int h = dof_offsets[i+1] - dof_offsets[i];
95 for (int j = 0; j < mat_r->NumColBlocks(); j++)
96 {
97 int w = dof_offsets[j+1] - dof_offsets[j];
98 mat_r->SetBlock(i,j,new SparseMatrix(h, w));
99 mat_i->SetBlock(i,j,new SparseMatrix(h, w));
100 }
101 }
102 y = new Vector(2*dof_offsets.Last());
103 *y=0.;
104 y_r = new BlockVector(*y, dof_offsets);
106}
107
109{
110 if (mat_r)
111 {
112 mat_r->Finalize(skip_zeros);
113 mat_i->Finalize(skip_zeros);
114 }
115 if (mat_e_r)
116 {
117 mat_e_r->Finalize(skip_zeros);
118 mat_e_i->Finalize(skip_zeros);
119 }
120 if (static_cond) { static_cond->Finalize(); }
121}
122
123/// Adds new Domain BF Integrator. Assumes ownership of @a bfi.
127 int n, int m)
128{
129 MFEM_VERIFY(n < trial_fes.Size(),
130 "ComplexDPGWeakFrom::AddTrialIntegrator: trial fespace index out of bounds");
131 MFEM_VERIFY(m < test_fecols.Size(),
132 "ComplexDPGWeakFrom::AddTrialIntegrator: test fecol index out of bounds");
133 if (bfi_r)
134 {
135 trial_integs_r(n,m)->Append(bfi_r);
136 }
137 if (bfi_i)
138 {
139 trial_integs_i(n,m)->Append(bfi_i);
140 }
141}
142
143/// Adds new Domain BF Integrator. Assumes ownership of @a bfi.
147 int n, int m)
148{
149 MFEM_VERIFY(n < test_fecols.Size() && m < test_fecols.Size(),
150 "ComplexDPGWeakFrom::AdTestIntegrator: test fecol index out of bounds");
151 if (bfi_r)
152 {
153 test_integs_r(n,m)->Append(bfi_r);
154 }
155 if (bfi_i)
156 {
157 test_integs_i(n,m)->Append(bfi_i);
158 }
159}
160
161/// Adds new Domain LF Integrator. Assumes ownership of @a bfi.
164 LinearFormIntegrator *lfi_i, int n)
165{
166 MFEM_VERIFY(n < test_fecols.Size(),
167 "ComplexDPGWeakFrom::AddDomainLFIntegrator: test fecol index out of bounds");
168 if (lfi_r)
169 {
170 lfis_r[n]->Append(lfi_r);
171 }
172 if (lfi_i)
173 {
174 lfis_i[n]->Append(lfi_i);
175 }
176}
177
179{
182 P->owns_blocks = 0;
183 R->owns_blocks = 0;
184 for (int i = 0; i<nblocks; i++)
185 {
186 const SparseMatrix *P_ = trial_fes[i]->GetConformingProlongation();
187 if (P_)
188 {
189 const SparseMatrix *R_ = trial_fes[i]->GetRestrictionMatrix();
190 P->SetBlock(i, i, const_cast<SparseMatrix*>(P_));
191 R->SetBlock(i, i, const_cast<SparseMatrix*>(R_));
192 }
193 }
194}
195
197{
198 Finalize(0);
199 if (!P) { BuildProlongation(); }
200
201 BlockMatrix * Pt = Transpose(*P);
202 BlockMatrix * PtA_r = mfem::Mult(*Pt, *mat_r);
203 BlockMatrix * PtA_i = mfem::Mult(*Pt, *mat_i);
204 mat_r->owns_blocks = 0;
205 mat_i->owns_blocks = 0;
206 for (int i = 0; i < nblocks; i++)
207 {
208 for (int j = 0; j < nblocks; j++)
209 {
210 SparseMatrix * tmp_r = &mat_r->GetBlock(i,j);
211 SparseMatrix * tmp_i = &mat_i->GetBlock(i,j);
212 if (Pt->IsZeroBlock(i, i))
213 {
214 PtA_r->SetBlock(i, j, tmp_r);
215 PtA_i->SetBlock(i, j, tmp_i);
216 }
217 else
218 {
219 delete tmp_r;
220 delete tmp_i;
221 }
222 }
223 }
224 delete mat_r;
225 delete mat_i;
226 if (mat_e_r)
227 {
228 BlockMatrix *PtAe_r = mfem::Mult(*Pt, *mat_e_r);
229 BlockMatrix *PtAe_i = mfem::Mult(*Pt, *mat_e_i);
230 mat_e_r->owns_blocks = 0;
231 mat_e_i->owns_blocks = 0;
232 for (int i = 0; i<nblocks; i++)
233 {
234 for (int j = 0; j<nblocks; j++)
235 {
236 SparseMatrix * tmp_r = &mat_e_r->GetBlock(i, j);
237 SparseMatrix * tmp_i = &mat_e_i->GetBlock(i, j);
238 if (Pt->IsZeroBlock(i, i))
239 {
240 PtAe_r->SetBlock(i, j, tmp_r);
241 PtAe_i->SetBlock(i, j, tmp_i);
242 }
243 else
244 {
245 delete tmp_r;
246 delete tmp_i;
247 }
248 }
249 }
250 delete mat_e_r;
251 delete mat_e_i;
252 mat_e_r = PtAe_r;
253 mat_e_i = PtAe_i;
254 }
255 delete Pt;
256
257 mat_r = mfem::Mult(*PtA_r, *P);
258 mat_i = mfem::Mult(*PtA_i, *P);
259
260 PtA_r->owns_blocks = 0;
261 PtA_i->owns_blocks = 0;
262 for (int i = 0; i < nblocks; i++)
263 {
264 for (int j = 0; j < nblocks; j++)
265 {
266 SparseMatrix * tmp_r = &PtA_r->GetBlock(j, i);
267 SparseMatrix * tmp_i = &PtA_i->GetBlock(j, i);
268 if (P->IsZeroBlock(i, i))
269 {
270 mat_r->SetBlock(j, i, tmp_r);
271 mat_i->SetBlock(j, i, tmp_i);
272 }
273 else
274 {
275 delete tmp_r;
276 delete tmp_i;
277 }
278 }
279 }
280 delete PtA_r;
281 delete PtA_i;
282
283 if (mat_e_r)
284 {
285 BlockMatrix *PtAeP_r = mfem::Mult(*mat_e_r, *P);
286 BlockMatrix *PtAeP_i = mfem::Mult(*mat_e_i, *P);
287 mat_e_r->owns_blocks = 0;
288 mat_e_i->owns_blocks = 0;
289 for (int i = 0; i < nblocks; i++)
290 {
291 for (int j = 0; j < nblocks; j++)
292 {
293 SparseMatrix * tmp_r = &mat_e_r->GetBlock(j, i);
294 SparseMatrix * tmp_i = &mat_e_i->GetBlock(j, i);
295 if (P->IsZeroBlock(i, i))
296 {
297 PtAeP_r->SetBlock(j, i, tmp_r);
298 PtAeP_i->SetBlock(j, i, tmp_i);
299 }
300 else
301 {
302 delete tmp_r;
303 delete tmp_i;
304 }
305 }
306 }
307
308 delete mat_e_r;
309 delete mat_e_i;
310 mat_e_r = PtAeP_r;
311 mat_e_i = PtAeP_i;
312 }
313 height = 2*mat_r->Height();
314 width = 2*mat_r->Width();
315}
316
317/// Assembles the form i.e. sums over all domain integrators.
319{
320 ElementTransformation *eltrans;
321 Array<int> faces, ori;
322
323 DofTransformation doftrans_i, doftrans_j;
324 if (mat_r == NULL)
325 {
326 AllocMat();
327 }
328
329 // loop through the elements
330 int dim = mesh->Dimension();
331 DenseMatrix B_r, Be_r, G_r, Ge_r, A_r;
332 DenseMatrix B_i, Be_i, G_i, Ge_i, A_i;
333 Vector vec_e_r, vec_r, b_r;
334 Vector vec_e_i, vec_i, b_i;
335 Array<int> vdofs;
336
337 // loop through elements
338 for (int iel = 0; iel < mesh -> GetNE(); iel++)
339 {
340 if (dim == 1)
341 {
342 mesh->GetElementVertices(iel, faces);
343 }
344 else if (dim == 2)
345 {
346 mesh->GetElementEdges(iel, faces, ori);
347 }
348 else if (dim == 3)
349 {
350 mesh->GetElementFaces(iel,faces,ori);
351 }
352 else
353 {
354 MFEM_ABORT("ComplexDPGWeakForm::Assemble: dim > 3 not supported");
355 }
356 int numfaces = faces.Size();
357
358 Array<int> test_offs(test_fecols.Size()+1); test_offs[0] = 0;
359 Array<int> trial_offs(trial_fes.Size()+1); trial_offs = 0;
360
361 eltrans = mesh->GetElementTransformation(iel);
362 for (int j = 0; j < test_fecols.Size(); j++)
363 {
364 int order = test_fecols[j]->GetOrder(); // assuming uniform order
365 test_offs[j+1] = test_fecols_vdims[j]*test_fecols[j]->GetFE(
366 eltrans->GetGeometryType(), order)->GetDof();
367 }
368 for (int j = 0; j < trial_fes.Size(); j++)
369 {
370 if (IsTraceFes[j])
371 {
372 for (int ie = 0; ie < faces.Size(); ie++)
373 {
374 trial_offs[j+1] += trial_fes[j]->GetVDim()*trial_fes[j]->GetFaceElement(
375 faces[ie])->GetDof();
376 }
377 }
378 else
379 {
380 trial_offs[j+1] = trial_fes[j]->GetVDim() * trial_fes[j]->GetFE(
381 iel)->GetDof();
382 }
383 }
384 test_offs.PartialSum();
385 trial_offs.PartialSum();
386
387 G_r.SetSize(test_offs.Last()); G_r = 0.0;
388 vec_r.SetSize(test_offs.Last()); vec_r = 0.0;
389 B_r.SetSize(test_offs.Last(),trial_offs.Last()); B_r = 0.0;
390 G_i.SetSize(test_offs.Last()); G_i = 0.0;
391 vec_i.SetSize(test_offs.Last()); vec_i = 0.0;
392 B_i.SetSize(test_offs.Last(),trial_offs.Last()); B_i = 0.0;
393
394 for (int j = 0; j < test_fecols.Size(); j++)
395 {
396 int order_j = test_fecols[j]->GetOrder();
397
398 eltrans = mesh->GetElementTransformation(iel);
399 const FiniteElement & test_fe =
400 *test_fecols[j]->GetFE(eltrans->GetGeometryType(), order_j);
401
402 // real integrators
403 for (int k = 0; k < lfis_r[j]->Size(); k++)
404 {
405 (*lfis_r[j])[k]->AssembleRHSElementVect(test_fe, *eltrans, vec_e_r);
406 vec_r.AddSubVector(vec_e_r, test_offs[j]);
407 }
408 // imag integrators
409 for (int k = 0; k < lfis_i[j]->Size(); k++)
410 {
411 (*lfis_i[j])[k]->AssembleRHSElementVect(test_fe,*eltrans,vec_e_i);
412 vec_i.AddSubVector(vec_e_i, test_offs[j]);
413 }
414
415 for (int i = 0; i < test_fecols.Size(); i++)
416 {
417 int order_i = test_fecols[i]->GetOrder();
418 eltrans = mesh->GetElementTransformation(iel);
419 const FiniteElement & test_fe_i =
420 *test_fecols[i]->GetFE(eltrans->GetGeometryType(), order_i);
421
422 // real integrators
423 for (int k = 0; k < test_integs_r(i,j)->Size(); k++)
424 {
425 if (i==j)
426 {
427 (*test_integs_r(i,j))[k]->AssembleElementMatrix(test_fe, *eltrans, Ge_r);
428 }
429 else
430 {
431 (*test_integs_r(i,j))[k]->AssembleElementMatrix2(test_fe_i, test_fe, *eltrans,
432 Ge_r);
433 }
434 G_r.AddSubMatrix(test_offs[j], test_offs[i], Ge_r);
435 }
436
437 // imag integrators
438 for (int k = 0; k < test_integs_i(i,j)->Size(); k++)
439 {
440 if (i==j)
441 {
442 (*test_integs_i(i,j))[k]->AssembleElementMatrix(test_fe,*eltrans,Ge_i);
443 }
444 else
445 {
446 (*test_integs_i(i,j))[k]->AssembleElementMatrix2(test_fe_i,test_fe,*eltrans,
447 Ge_i);
448 }
449 G_i.AddSubMatrix(test_offs[j], test_offs[i], Ge_i);
450 }
451 }
452
453 for (int i = 0; i < trial_fes.Size(); i++)
454 {
455 if (IsTraceFes[i])
456 {
457 // real integrators
458 for (int k = 0; k < trial_integs_r(i,j)->Size(); k++)
459 {
460 int face_dof_offs = 0;
461 for (int ie = 0; ie < numfaces; ie++)
462 {
463 int iface = faces[ie];
465 const FiniteElement & tfe = *trial_fes[i]->GetFaceElement(iface);
466 (*trial_integs_r(i,j))[k]->AssembleTraceFaceMatrix(iel, tfe, test_fe, *ftr,
467 Be_r);
468 B_r.AddSubMatrix(test_offs[j], trial_offs[i]+face_dof_offs, Be_r);
469 face_dof_offs += Be_r.Width();
470 }
471 }
472 // imag integrators
473 for (int k = 0; k < trial_integs_i(i,j)->Size(); k++)
474 {
475 int face_dof_offs = 0;
476 for (int ie = 0; ie < numfaces; ie++)
477 {
478 int iface = faces[ie];
480 const FiniteElement & tfe = *trial_fes[i]->GetFaceElement(iface);
481 (*trial_integs_i(i,j))[k]->AssembleTraceFaceMatrix(iel,tfe,test_fe,*ftr,Be_i);
482 B_i.AddSubMatrix(test_offs[j], trial_offs[i]+face_dof_offs, Be_i);
483 face_dof_offs += Be_i.Width();
484 }
485 }
486 }
487 else
488 {
489 const FiniteElement & fe = *trial_fes[i]->GetFE(iel);
490 eltrans = mesh->GetElementTransformation(iel);
491 // real integrators
492 for (int k = 0; k < trial_integs_r(i,j)->Size(); k++)
493 {
494 (*trial_integs_r(i,j))[k]->AssembleElementMatrix2(fe,test_fe,*eltrans,Be_r);
495 B_r.AddSubMatrix(test_offs[j], trial_offs[i], Be_r);
496 }
497 // imag integrators
498 for (int k = 0; k < trial_integs_i(i,j)->Size(); k++)
499 {
500 (*trial_integs_i(i,j))[k]->AssembleElementMatrix2(fe, test_fe, *eltrans, Be_i);
501 B_i.AddSubMatrix(test_offs[j], trial_offs[i], Be_i);
502 }
503 }
504 }
505 }
506
507 ComplexCholeskyFactors chol(G_r.GetData(), G_i.GetData());
508 int h = G_r.Height();
509 chol.Factor(h);
510
511 int w = B_r.Width();
512 chol.LSolve(h,w,B_r.GetData(), B_i.GetData());
513 chol.LSolve(h,1,vec_r.GetData(), vec_i.GetData());
514
515 Vector vec(vec_i.Size()+vec_r.Size());
516 vec.SetVector(vec_r, 0);
517 vec.SetVector(vec_i, vec_r.Size());
518
519 if (store_matrices)
520 {
521 Bmat[iel] = new ComplexDenseMatrix(new DenseMatrix(B_r), new DenseMatrix(B_i),
522 true,true);
523 fvec[iel] = new Vector(vec);
524 }
525 ComplexDenseMatrix B(&B_r, &B_i, false, false);
527 Vector b(B.Width());
528 B.MultTranspose(vec, b);
529
530 b_r.MakeRef(b, 0, b.Size()/2);
531 b_i.MakeRef(b, b.Size()/2,b.Size()/2);
532
533 if (static_cond)
534 {
535 static_cond->AssembleReducedSystem(iel,*A,b_r,b_i);
536 }
537 else
538 {
539 // Assembly
540 for (int i = 0; i<trial_fes.Size(); i++)
541 {
542 Array<int> vdofs_i;
543 doftrans_i.SetDofTransformation(nullptr);
544 if (IsTraceFes[i])
545 {
546 Array<int> face_vdofs;
547 for (int k = 0; k < numfaces; k++)
548 {
549 int iface = faces[k];
550 trial_fes[i]->GetFaceVDofs(iface, face_vdofs);
551 vdofs_i.Append(face_vdofs);
552 }
553 }
554 else
555 {
556 trial_fes[i]->GetElementVDofs(iel, vdofs_i, doftrans_i);
557 }
558 for (int j = 0; j < trial_fes.Size(); j++)
559 {
560 Array<int> vdofs_j;
561 doftrans_j.SetDofTransformation(nullptr);
562
563 if (IsTraceFes[j])
564 {
565 Array<int> face_vdofs;
566 for (int k = 0; k < numfaces; k++)
567 {
568 int iface = faces[k];
569 trial_fes[j]->GetFaceVDofs(iface, face_vdofs);
570 vdofs_j.Append(face_vdofs);
571 }
572 }
573 else
574 {
575 trial_fes[j]->GetElementVDofs(iel, vdofs_j, doftrans_j);
576 }
577
578 DenseMatrix Ae_r, Ae_i;
579 A->real().GetSubMatrix(trial_offs[i],trial_offs[i+1],
580 trial_offs[j],trial_offs[j+1], Ae_r);
581 A->imag().GetSubMatrix(trial_offs[i],trial_offs[i+1],
582 trial_offs[j],trial_offs[j+1], Ae_i);
583 TransformDual(doftrans_i, doftrans_j, Ae_r);
584 TransformDual(doftrans_i, doftrans_j, Ae_i);
585 if (!mat_r)
586 {
587 mfem::out << "null matrix " << std::endl;
588 }
589 mat_r->GetBlock(i,j).AddSubMatrix(vdofs_i,vdofs_j, Ae_r);
590 mat_i->GetBlock(i,j).AddSubMatrix(vdofs_i,vdofs_j, Ae_i);
591 }
592
593 // assemble rhs
594 Vector vec1_r(b_r,trial_offs[i],trial_offs[i+1]-trial_offs[i]);
595 Vector vec1_i(b_i,trial_offs[i],trial_offs[i+1]-trial_offs[i]);
596
597 doftrans_i.TransformDual(vec1_r);
598 doftrans_i.TransformDual(vec1_i);
599
600 y_r->GetBlock(i).AddElementVector(vdofs_i,vec1_r);
601 y_i->GetBlock(i).AddElementVector(vdofs_i,vec1_i);
602 }
603 }
604 delete A;
605 } // end of loop through elements
606}
607
609 &ess_tdof_list,
610 Vector &x,
612 Vector &X,
613 Vector &B,
614 int copy_interior)
615{
616 FormSystemMatrix(ess_tdof_list, A);
617 if (static_cond)
618 {
619 static_cond->ReduceSystem(x, X, B, copy_interior);
620 }
621 else if (!P)
622 {
623 Vector x_r(x, 0, x.Size()/2);
624 Vector x_i(x, x.Size()/2, x.Size()/2);
625 EliminateVDofsInRHS(ess_tdof_list, x_r,x_i, *y_r, *y_i);
626 if (!copy_interior)
627 {
628 x_r.SetSubVectorComplement(ess_tdof_list, 0.0);
629 x_i.SetSubVectorComplement(ess_tdof_list, 0.0);
630 }
631 X.MakeRef(x, 0, x.Size());
632 B.MakeRef(*y,0,y->Size());
633 }
634 else // non conforming space
635 {
636 B.SetSize(2*P->Width());
637 Vector B_r(B, 0, P->Width());
638 Vector B_i(B, P->Width(),P->Width());
639
640 P->MultTranspose(*y_r, B_r);
641 P->MultTranspose(*y_i, B_i);
642 Vector tmp_r,tmp_i;
643 for (int i = 0; i<nblocks; i++)
644 {
645 if (P->IsZeroBlock(i,i))
646 {
647 int offset = tdof_offsets[i];
648 tmp_r.MakeRef(*y_r, offset,tdof_offsets[i+1]-tdof_offsets[i]);
649 tmp_i.MakeRef(*y_i, offset,tdof_offsets[i+1]-tdof_offsets[i]);
650 B_r.SetVector(tmp_r,offset);
651 B_i.SetVector(tmp_i,offset);
652 }
653 }
654
655 X.SetSize(2*R->Height());
656 Vector X_r(X, 0, X.Size()/2);
657 Vector X_i(X, X.Size()/2, X.Size()/2);
658
659 Vector x_r(x, 0,x.Size()/2);
660 Vector x_i(x, x.Size()/2, x.Size()/2);
661
662 R->Mult(x_r, X_r);
663 R->Mult(x_i, X_i);
664 for (int i = 0; i<nblocks; i++)
665 {
666 if (R->IsZeroBlock(i,i))
667 {
668 int offset = tdof_offsets[i];
669 tmp_r.MakeRef(x_r, offset, tdof_offsets[i+1]-tdof_offsets[i]);
670 tmp_i.MakeRef(x_i, offset, tdof_offsets[i+1]-tdof_offsets[i]);
671 X_r.SetVector(tmp_r,offset);
672 X_i.SetVector(tmp_i,offset);
673 }
674 }
675
676 EliminateVDofsInRHS(ess_tdof_list, X_r, X_i, B_r, B_i);
677 if (!copy_interior)
678 {
679 X_r.SetSubVectorComplement(ess_tdof_list, 0.0);
680 X_i.SetSubVectorComplement(ess_tdof_list, 0.0);
681 }
682 }
683}
684
686 &ess_tdof_list,
688{
689 if (static_cond)
690 {
692 {
693 static_cond->SetEssentialTrueDofs(ess_tdof_list);
695 }
697 }
698 else
699 {
700 if (!mat_e_r)
701 {
702 bool conforming = true;
703 for (int i = 0; i<nblocks; i++)
704 {
705 const SparseMatrix *P_ = trial_fes[i]->GetConformingProlongation();
706 if (P_)
707 {
708 conforming = false;
709 break;
710 }
711 }
712 if (!conforming) { ConformingAssemble(); }
713 const int remove_zeros = 0;
714 EliminateVDofs(ess_tdof_list, diag_policy);
715 Finalize(remove_zeros);
716 }
717 mat = new ComplexOperator(mat_r,mat_i,false,false);
718 A.Reset(mat,false);
719 }
720}
721
723 const Array<int> &vdofs, const Vector &x_r, const Vector & x_i,
724 Vector &b_r, Vector & b_i)
725{
726 mat_e_r->AddMult(x_r,b_r,-1.);
727 mat_e_i->AddMult(x_i,b_r,1.);
728 mat_e_r->AddMult(x_i,b_i,-1.);
729 mat_e_i->AddMult(x_r,b_i,-1.);
730
731 mat_r->PartMult(vdofs,x_r,b_r);
732 mat_r->PartMult(vdofs,x_i,b_i);
733}
734
737{
738 if (mat_e_r == NULL)
739 {
740 Array<int> offsets;
741
742 offsets.MakeRef( (P) ? tdof_offsets : dof_offsets);
743
744 mat_e_r = new BlockMatrix(offsets);
745 mat_e_r->owns_blocks = 1;
746 mat_e_i = new BlockMatrix(offsets);
747 mat_e_i->owns_blocks = 1;
748 for (int i = 0; i < mat_e_r->NumRowBlocks(); i++)
749 {
750 int h = offsets[i+1] - offsets[i];
751 for (int j = 0; j < mat_e_r->NumColBlocks(); j++)
752 {
753 int w = offsets[j+1] - offsets[j];
754 mat_e_r->SetBlock(i, j, new SparseMatrix(h, w));
755 mat_e_i->SetBlock(i, j, new SparseMatrix(h, w));
756 }
757 }
758 }
761}
762
764 Vector &x)
765{
766 if (static_cond)
767 {
768 // Private dofs back solve
770 }
771 else if (!P)
772 {
773 x.SyncMemory(X);
774 }
775 else
776 {
777 x.SetSize(2*P->Height());
778 Vector X_r(const_cast<Vector &>(X), 0, X.Size()/2);
779 Vector X_i(const_cast<Vector &>(X), X.Size()/2, X.Size()/2);
780
781 Vector x_r(x, 0, x.Size()/2);
782 Vector x_i(x, x.Size()/2, x.Size()/2);
783
784 P->Mult(X_r, x_r);
785 P->Mult(X_i, x_i);
786
787 Vector tmp_r, tmp_i;
788 for (int i = 0; i<nblocks; i++)
789 {
790 if (P->IsZeroBlock(i,i))
791 {
792 int offset = tdof_offsets[i];
793 tmp_r.MakeRef(X_r, offset, tdof_offsets[i+1]-tdof_offsets[i]);
794 tmp_i.MakeRef(X_i, offset, tdof_offsets[i+1]-tdof_offsets[i]);
795 x_r.SetVector(tmp_r,offset);
796 x_i.SetVector(tmp_i,offset);
797 }
798 }
799 }
800}
801
803{
804 if (initialized)
805 {
806 for (int k = 0; k < trial_integs_r.NumRows(); k++)
807 {
808 for (int l = 0; l < trial_integs_r.NumCols(); l++)
809 {
810 for (int i = 0; i < trial_integs_r(k,l)->Size(); i++)
811 {
812 delete (*trial_integs_r(k,l))[i];
813 }
814 delete trial_integs_r(k,l);
815 for (int i = 0; i < trial_integs_i(k,l)->Size(); i++)
816 {
817 delete (*trial_integs_i(k,l))[i];
818 }
819 delete trial_integs_i(k,l);
820 }
821 }
822 trial_integs_r.DeleteAll();
823 trial_integs_i.DeleteAll();
824
825 for (int k = 0; k < test_integs_r.NumRows(); k++)
826 {
827 for (int l = 0; l < test_integs_r.NumCols(); l++)
828 {
829 for (int i = 0; i < test_integs_r(k,l)->Size(); i++)
830 {
831 delete (*test_integs_r(k,l))[i];
832 }
833 delete test_integs_r(k,l);
834 for (int i = 0; i < test_integs_i(k,l)->Size(); i++)
835 {
836 delete (*test_integs_i(k,l))[i];
837 }
838 delete test_integs_i(k,l);
839 }
840 }
841 test_integs_r.DeleteAll();
842 test_integs_i.DeleteAll();
843
844 for (int k = 0; k < lfis_r.Size(); k++)
845 {
846 for (int i = 0; i < lfis_r[k]->Size(); i++)
847 {
848 delete (*lfis_r[k])[i];
849 }
850 delete lfis_r[k];
851 for (int i = 0; i < lfis_i[k]->Size(); i++)
852 {
853 delete (*lfis_i[k])[i];
854 }
855 delete lfis_i[k];
856 }
857 lfis_r.DeleteAll();
858 lfis_i.DeleteAll();
859 }
860}
861
863{
864 delete mat_e_r; mat_e_r = nullptr;
865 delete mat_e_i; mat_e_i = nullptr;
866 delete mat; mat = nullptr;
867 delete mat_r; mat_r = nullptr;
868 delete mat_i; mat_i = nullptr;
869 delete y; y = nullptr;
870 delete y_r; y_r = nullptr;
871 delete y_i; y_i = nullptr;
872
873 if (P)
874 {
875 delete P; P = nullptr;
876 delete R; R = nullptr;
877 }
878
879 if (static_cond)
880 {
882 }
883 else
884 {
885 delete static_cond; static_cond = nullptr;
886 }
887
889
892 width = height;
893
894 initialized = true;
895
896 if (store_matrices)
897 {
898 for (int i = 0; i < Bmat.Size(); i++)
899 {
900 delete Bmat[i]; Bmat[i] = nullptr;
901 delete fvec[i]; fvec[i] = nullptr;
902 }
904 fvec.SetSize(mesh->GetNE());
905 for (int i = 0; i < Bmat.Size(); i++)
906 {
907 Bmat[i] = nullptr;
908 fvec[i] = nullptr;
909 }
910 }
911}
912
918
920{
921 MFEM_VERIFY(store_matrices,
922 "Matrices needed for the residual are not store. Call ComplexDPGWeakForm::StoreMatrices()")
923 // wrap vector in a blockvector
924 int n = x.Size()/2;
925
926 BlockVector x_r(const_cast<Vector &>(x),0,dof_offsets);
927 BlockVector x_i(const_cast<Vector &>(x),n,dof_offsets);
928
929 // Element vector of trial space size
930 Vector u;
931 Array<int> vdofs;
932 Array<int> faces, ori;
933 int dim = mesh->Dimension();
935 // loop through elements
936 for (int iel = 0; iel < mesh -> GetNE(); iel++)
937 {
938 if (dim == 1)
939 {
940 mesh->GetElementVertices(iel, faces);
941 }
942 else if (dim == 2)
943 {
944 mesh->GetElementEdges(iel, faces, ori);
945 }
946 else if (dim == 3)
947 {
948 mesh->GetElementFaces(iel,faces,ori);
949 }
950 else
951 {
952 MFEM_ABORT("ComplexDPGWeakForm::ComputeResidual: "
953 "dim > 3 not supported");
954 }
955 int numfaces = faces.Size();
956
957 Array<int> trial_offs(trial_fes.Size()+1); trial_offs = 0;
958
959 for (int j = 0; j < trial_fes.Size(); j++)
960 {
961 if (IsTraceFes[j])
962 {
963 for (int ie = 0; ie < faces.Size(); ie++)
964 {
965 trial_offs[j+1] += trial_fes[j]->GetFaceElement(faces[ie])->GetDof();
966 }
967 }
968 else
969 {
970 trial_offs[j+1] = trial_fes[j]->GetVDim() * trial_fes[j]->GetFE(
971 iel)->GetDof();
972 }
973 }
974 trial_offs.PartialSum();
975
976 int nn = trial_offs.Last();
977 u.SetSize(2*nn);
978 DofTransformation doftrans;
979 for (int i = 0; i<trial_fes.Size(); i++)
980 {
981 vdofs.SetSize(0);
982 doftrans.SetDofTransformation(nullptr);
983 if (IsTraceFes[i])
984 {
985 Array<int> face_vdofs;
986 for (int k = 0; k < numfaces; k++)
987 {
988 int iface = faces[k];
989 trial_fes[i]->GetFaceVDofs(iface, face_vdofs);
990 vdofs.Append(face_vdofs);
991 }
992 }
993 else
994 {
995 trial_fes[i]->GetElementVDofs(iel, vdofs, doftrans);
996 }
997 Vector vec1_r;
998 Vector vec1_i;
999 vec1_r.MakeRef(u, trial_offs[i], trial_offs[i+1]-trial_offs[i]);
1000 vec1_i.MakeRef(u, trial_offs[i]+nn, trial_offs[i+1]-trial_offs[i]);
1001 x_r.GetBlock(i).GetSubVector(vdofs,vec1_r);
1002 x_i.GetBlock(i).GetSubVector(vdofs,vec1_i);
1003 doftrans.InvTransformPrimal(vec1_r);
1004 doftrans.InvTransformPrimal(vec1_i);
1005 } // end of loop through trial spaces
1006
1007 // residual
1008 Vector v(Bmat[iel]->Height());
1009 Bmat[iel]->Mult(u,v);
1010 v -= *fvec[iel];
1011 residuals[iel] = v.Norml2();
1012 } // end of loop through elements
1013 return residuals;
1014}
1015
1017{
1018 delete mat_e_r; mat_e_r = nullptr;
1019 delete mat_e_i; mat_e_i = nullptr;
1020 delete mat; mat = nullptr;
1021 delete mat_r; mat_r = nullptr;
1022 delete mat_i; mat_i = nullptr;
1023 delete y; y = nullptr;
1024 delete y_r; y_r = nullptr;
1025 delete y_i; y_i = nullptr;
1026
1028
1029 if (P)
1030 {
1031 delete P;
1032 delete R;
1033 }
1034
1035 delete static_cond;
1036
1037 if (store_matrices)
1038 {
1039 for (int i = 0; i<mesh->GetNE(); i++)
1040 {
1041 delete Bmat[i]; Bmat[i] = nullptr;
1042 delete fvec[i]; fvec[i] = nullptr;
1043 }
1044 }
1045}
1046
1047} // namespace mfem
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:840
int Size() const
Return the logical size of the array.
Definition array.hpp:166
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
Definition array.cpp:104
void MakeRef(T *data_, int size_, bool own_data=false)
Make this Array a reference to a pointer.
Definition array.hpp:1053
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:912
T & Last()
Return the last element in the array.
Definition array.hpp:945
Abstract base class BilinearFormIntegrator.
void MultTranspose(const Vector &x, Vector &y) const override
MatrixTranspose-Vector Multiplication y = A'*x.
void PartMult(const Array< int > &rows, const Vector &x, Vector &y) const
Partial matrix vector multiplication of (*this) with x involving only the rows given by rows....
void SetBlock(int i, int j, SparseMatrix *mat)
Set A(i,j) = mat.
void Mult(const Vector &x, Vector &y) const override
Matrix-Vector Multiplication y = A*x.
int NumColBlocks() const
Return the number of column blocks.
int IsZeroBlock(int i, int j) const
Check if block (i,j) is a zero block.
void EliminateRowCols(const Array< int > &vdofs, BlockMatrix *Ae, DiagonalPolicy dpolicy=DIAG_ONE)
Eliminate the rows and columns corresponding to the entries in vdofs + save the eliminated entries in...
SparseMatrix & GetBlock(int i, int j)
Return a reference to block (i,j). Reference may be invalid if Aij(i,j) == NULL.
void AddMult(const Vector &x, Vector &y, const real_t val=1.) const override
Matrix-Vector Multiplication y = y + val*A*x.
int owns_blocks
If owns_blocks the SparseMatrix objects Aij will be deallocated.
void Finalize(int skip_zeros=1) override
Finalize all the submatrices.
int NumRowBlocks() const
Return the number of row blocks.
A class to handle Vectors in a block fashion.
Vector & GetBlock(int i)
Get the i-th vector in the block.
Class that performs static condensation of interior dofs for multiple FE spaces for complex systems (...
void Finalize(int skip_zeros=0)
Finalize the construction of the Schur complement matrix.
void AssembleReducedSystem(int el, ComplexDenseMatrix &elmat, Vector &elvect_r, Vector &elvect_i)
void FormSystemMatrix(Operator::DiagonalPolicy diag_policy)
void SetEssentialTrueDofs(const Array< int > &ess_tdof_list)
Determine and save internally essential reduced true dofs.
void ComputeSolution(const Vector &sc_sol, Vector &sol) const
void ReduceSystem(Vector &x, Vector &X, Vector &B, int copy_interior=0) const
Set the reduced solution X and r.h.s B vectors from the full linear system solution x and r....
void LSolve(int m, int n, real_t *X_r, real_t *X_i) const
bool Factor(int m, real_t TOL=0.0) override
Compute the Cholesky factorization of the current matrix.
void AddTestIntegrator(BilinearFormIntegrator *bfi_r, BilinearFormIntegrator *bfi_i, int n, int m)
Adds new Test Integrator. Assumes ownership of bfi_r and bfi_i.
BlockMatrix * mat_e_r
Block Matrix used to store the eliminations from the b.c. Owned. .
Array2D< Array< BilinearFormIntegrator * > * > test_integs_r
Set of Test Space (broken) Integrators to be applied for matrix G.
Array< FiniteElementSpace * > trial_fes
Trial FE spaces.
Array< int > IsTraceFes
Flags to determine if a FiniteElementSpace is Trace.
Array< ComplexDenseMatrix * > Bmat
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Array2D< Array< BilinearFormIntegrator * > * > test_integs_i
mfem::Operator::DiagonalPolicy diag_policy
void AddTrialIntegrator(BilinearFormIntegrator *bfi_r, BilinearFormIntegrator *bfi_i, int n, int m)
Adds new Domain BF Integrator. Assumes ownership of bfi.
void EliminateVDofsInRHS(const Array< int > &vdofs, const Vector &x_r, const Vector &x_i, Vector &b_r, Vector &b_i)
BlockMatrix * P
Block Prolongation.
BlockMatrix * R
Block Restriction.
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Array2D< Array< BilinearFormIntegrator * > * > trial_integs_r
Set of Trial Integrators to be applied for matrix B.
void AddDomainLFIntegrator(LinearFormIntegrator *lfi_r, LinearFormIntegrator *lfi_i, int n)
Adds new Domain LF Integrator. Assumes ownership of lfi_r and lfi_i.
BlockMatrix * mat_r
Block matrix to be associated with the real/imag Block bilinear form. Owned.
Array2D< Array< BilinearFormIntegrator * > * > trial_integs_i
void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
Array< Array< LinearFormIntegrator * > * > lfis_i
Vector & ComputeResidual(const Vector &x)
void EliminateVDofs(const Array< int > &vdofs, Operator::DiagonalPolicy dpolicy=Operator::DIAG_ONE)
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all integrators.
virtual ~ComplexDPGWeakForm()
Destroys bilinear form.
Array< FiniteElementCollection * > test_fecols
Test FE Collections (Broken)
BlockVector * y_r
BlockVectors to be associated with the real/imag Block linear form.
Array< Array< LinearFormIntegrator * > * > lfis_r
Set of LinearForm Integrators to be applied.
virtual void RecoverFEMSolution(const Vector &X, Vector &x)
ComplexBlockStaticCondensation * static_cond
Owned.
Specialization of the ComplexOperator built from a pair of Dense Matrices. The purpose of this specia...
DenseMatrix & imag() override
DenseMatrix & real() override
Real or imaginary part accessor methods.
Mimic the action of a complex operator using two real operators.
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void GetSubMatrix(const Array< int > &idx, DenseMatrix &A) const
real_t * GetData() const
Returns the matrix data array. Warning: this method casts away constness.
Definition densemat.hpp:126
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:116
void AddSubMatrix(const Array< int > &idx, const DenseMatrix &A)
(*this)(idx[i],idx[j]) += A(i,j)
void TransformDual(real_t *v) const
Definition doftrans.cpp:77
void InvTransformPrimal(real_t *v) const
Definition doftrans.cpp:47
void SetDofTransformation(const StatelessDofTransformation &dof_trans)
Set or change the nested StatelessDofTransformation object.
Definition doftrans.hpp:176
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition eltrans.hpp:175
A specialized ElementTransformation class representing a face and its two neighboring elements.
Definition eltrans.hpp:750
Abstract class for all finite elements.
Definition fe_base.hpp:247
Abstract base class LinearFormIntegrator.
Definition lininteg.hpp:28
virtual FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
Definition mesh.cpp:1104
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition mesh.hpp:1609
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1377
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1306
void GetElementTransformation(int i, IsoparametricTransformation *ElTr) const
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and wi...
Definition mesh.cpp:360
void GetElementFaces(int i, Array< int > &faces, Array< int > &ori) const
Return the indices and the orientations of all faces of element i.
Definition mesh.cpp:7835
void GetElementEdges(int i, Array< int > &edges, Array< int > &cor) const
Return the indices and the orientations of all edges of element i.
Definition mesh.cpp:7588
Pointer to an Operator of a specified type.
Definition handle.hpp:34
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition handle.hpp:145
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
DiagonalPolicy
Defines operator diagonal policy upon elimination of rows and/or columns.
Definition operator.hpp:48
@ DIAG_ONE
Set the diagonal value to one.
Definition operator.hpp:50
@ DIAG_ZERO
Set the diagonal value to zero.
Definition operator.hpp:49
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition operator.hpp:72
Data type sparse matrix.
Definition sparsemat.hpp:51
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Vector data type.
Definition vector.hpp:82
void SetVector(const Vector &v, int offset)
(*this)[i + offset] = v[i]
Definition vector.cpp:353
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add elements of the elemvect Vector to the entries listed in dofs. Negative dof values cause the -dof...
Definition vector.cpp:785
void AddSubVector(const Vector &v, int offset)
(*this)[i + offset] += v[i]
Definition vector.cpp:365
real_t Norml2() const
Returns the l2 norm of the vector.
Definition vector.cpp:968
void SyncMemory(const Vector &v) const
Update the memory location of the vector to match v.
Definition vector.hpp:272
int Size() const
Returns the size of the vector.
Definition vector.hpp:234
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:584
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
Definition vector.hpp:243
void SetSubVectorComplement(const Array< int > &dofs, const real_t val)
Set all vector entries NOT in the dofs Array to the given val.
Definition vector.cpp:854
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition vector.cpp:676
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition vector.hpp:660
int dim
Definition ex24.cpp:53
real_t b
Definition lissajous.cpp:42
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition table.cpp:505
void TransformDual(const DofTransformation &ran_dof_trans, const DofTransformation &dom_dof_trans, DenseMatrix &elmat)
Definition doftrans.cpp:152
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition globals.hpp:66
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
Definition table.cpp:443
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.