MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
nonlinearform.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 "fem.hpp"
13#include "../general/forall.hpp"
14
15namespace mfem
16{
17
19{
20 if (ext)
21 {
22 MFEM_ABORT("the assembly level has already been set!");
23 }
24 assembly = assembly_level;
25 switch (assembly)
26 {
28 ext = new MFNonlinearFormExtension(this);
29 break;
31 ext = new PANonlinearFormExtension(this);
32 break;
34 // This is the default
35 break;
36 default:
37 mfem_error("Unknown assembly level for this form.");
38 }
39}
40
41void NonlinearForm::SetEssentialBC(const Array<int> &bdr_attr_is_ess,
42 Vector *rhs)
43{
44 // virtual call, works in parallel too
45 fes->GetEssentialTrueDofs(bdr_attr_is_ess, ess_tdof_list);
46
47 if (rhs)
48 {
49 for (int i = 0; i < ess_tdof_list.Size(); i++)
50 {
51 (*rhs)(ess_tdof_list[i]) = 0.0;
52 }
53 }
54}
55
57{
58 if (!P)
59 {
60 ess_vdofs_list.Copy(ess_tdof_list); // ess_vdofs_list --> ess_tdof_list
61 }
62 else
63 {
64 Array<int> ess_vdof_marker, ess_tdof_marker;
66 ess_vdof_marker);
67 if (Serial())
68 {
69 fes->ConvertToConformingVDofs(ess_vdof_marker, ess_tdof_marker);
70 }
71 else
72 {
73#ifdef MFEM_USE_MPI
75 ess_tdof_marker.SetSize(pf->GetTrueVSize());
76 pf->Dof_TrueDof_Matrix()->BooleanMultTranspose(1, ess_vdof_marker,
77 0, ess_tdof_marker);
78#else
79 MFEM_ABORT("internal MFEM error");
80#endif
81 }
83 }
84}
85
87{
88 if (ext)
89 {
90 MFEM_VERIFY(!fnfi.Size(), "Interior faces terms not yet implemented!");
91 MFEM_VERIFY(!bfnfi.Size(), "Boundary face terms not yet implemented!");
92 return ext->GetGridFunctionEnergy(x);
93 }
94
95 Array<int> vdofs;
96 Vector el_x;
97 const FiniteElement *fe;
99 Mesh *mesh = fes->GetMesh();
100 real_t energy = 0.0;
101
102 if (dnfi.Size())
103 {
104 // Which attributes need to be processed?
105 Array<int> attr_marker(mesh->attributes.Size() ?
106 mesh->attributes.Max() : 0);
107 attr_marker = 0;
108 for (int k = 0; k < dnfi.Size(); k++)
109 {
110 if (dnfi_marker[k] == NULL)
111 {
112 attr_marker = 1;
113 break;
114 }
115 Array<int> &marker = *dnfi_marker[k];
116 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
117 "invalid marker for domain integrator #"
118 << k << ", counting from zero");
119 for (int i = 0; i < attr_marker.Size(); i++)
120 {
121 attr_marker[i] |= marker[i];
122 }
123 }
124
125 DofTransformation doftrans;
126 for (int i = 0; i < fes->GetNE(); i++)
127 {
128 const int attr = mesh->GetAttribute(i);
129 if (attr_marker[attr-1] == 0) { continue; }
130
131 fe = fes->GetFE(i);
132 fes->GetElementVDofs(i, vdofs, doftrans);
134 x.GetSubVector(vdofs, el_x);
135 doftrans.InvTransformPrimal(el_x);
136 for (int k = 0; k < dnfi.Size(); k++)
137 {
138 if (dnfi_marker[k] &&
139 (*dnfi_marker[k])[attr-1] == 0) { continue; }
140
141 energy += dnfi[k]->GetElementEnergy(*fe, *T, el_x);
142 }
143 }
144 }
145
146 if (bnfi.Size())
147 {
148 // Which boundary attributes need to be processed?
149 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
150 mesh->bdr_attributes.Max() : 0);
151 bdr_attr_marker = 0;
152 for (int k = 0; k < bnfi.Size(); k++)
153 {
154 if (bnfi_marker[k] == NULL)
155 {
156 bdr_attr_marker = 1;
157 break;
158 }
159 Array<int> &bdr_marker = *bnfi_marker[k];
160 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
161 "invalid boundary marker for boundary integrator #"
162 << k << ", counting from zero");
163 for (int i = 0; i < bdr_attr_marker.Size(); i++)
164 {
165 bdr_attr_marker[i] |= bdr_marker[i];
166 }
167 }
168
169 DofTransformation doftrans;
170 for (int i = 0; i < fes->GetNBE(); i++)
171 {
172 const int bdr_attr = mesh->GetBdrAttribute(i);
173 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
174
175 fe = fes->GetBE(i);
176 fes->GetBdrElementVDofs(i, vdofs, doftrans);
178 x.GetSubVector(vdofs, el_x);
179 doftrans.InvTransformPrimal(el_x);
180 for (int k = 0; k < bnfi.Size(); k++)
181 {
182 if (bnfi_marker[k] &&
183 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
184
185 energy += bnfi[k]->GetElementEnergy(*fe, *T, el_x);
186 }
187 }
188
189 }
190
191 if (fnfi.Size())
192 {
193 MFEM_ABORT("TODO: add energy contribution from interior face terms");
194 }
195
196 if (bfnfi.Size())
197 {
198 MFEM_ABORT("TODO: add energy contribution from boundary face terms");
199 }
200
201 return energy;
202}
203
205{
206 MFEM_VERIFY(x.Size() == Width(), "invalid input Vector size");
207 if (P)
208 {
209 aux1.SetSize(P->Height());
210 P->Mult(x, aux1);
211 return aux1;
212 }
213 return x;
214}
215
216void NonlinearForm::Mult(const Vector &x, Vector &y) const
217{
218 const Vector &px = Prolongate(x);
219 if (P) { aux2.SetSize(P->Height()); }
220
221 // If we are in parallel, ParNonLinearForm::Mult uses the aux2 vector. In
222 // serial, place the result directly in y (when there is no P).
223 Vector &py = P ? aux2 : y;
224
225 if (ext)
226 {
227 ext->Mult(px, py);
228 if (Serial())
229 {
230 if (cP) { cP->MultTranspose(py, y); }
231 const int N = ess_tdof_list.Size();
232 const auto tdof = ess_tdof_list.Read();
233 auto Y = y.ReadWrite();
234 mfem::forall(N, [=] MFEM_HOST_DEVICE (int i) { Y[tdof[i]] = 0.0; });
235 }
236 // In parallel, the result is in 'py' which is an alias for 'aux2'.
237 return;
238 }
239
240 Array<int> vdofs;
241 Vector el_x, el_y;
242 const FiniteElement *fe;
244 Mesh *mesh = fes->GetMesh();
245
246 py = 0.0;
247
248 if (dnfi.Size())
249 {
250 // Which attributes need to be processed?
251 Array<int> attr_marker(mesh->attributes.Size() ?
252 mesh->attributes.Max() : 0);
253 attr_marker = 0;
254 for (int k = 0; k < dnfi.Size(); k++)
255 {
256 if (dnfi_marker[k] == NULL)
257 {
258 attr_marker = 1;
259 break;
260 }
261 Array<int> &marker = *dnfi_marker[k];
262 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
263 "invalid marker for domain integrator #"
264 << k << ", counting from zero");
265 for (int i = 0; i < attr_marker.Size(); i++)
266 {
267 attr_marker[i] |= marker[i];
268 }
269 }
270
271 DofTransformation doftrans;
272 for (int i = 0; i < fes->GetNE(); i++)
273 {
274 const int attr = mesh->GetAttribute(i);
275 if (attr_marker[attr-1] == 0) { continue; }
276
277 fe = fes->GetFE(i);
278 fes->GetElementVDofs(i, vdofs, doftrans);
280 px.GetSubVector(vdofs, el_x);
281 doftrans.InvTransformPrimal(el_x);
282 for (int k = 0; k < dnfi.Size(); k++)
283 {
284 if (dnfi_marker[k] &&
285 (*dnfi_marker[k])[attr-1] == 0) { continue; }
286
287 dnfi[k]->AssembleElementVector(*fe, *T, el_x, el_y);
288 doftrans.TransformDual(el_y);
289 py.AddElementVector(vdofs, el_y);
290 }
291 }
292 }
293
294 if (bnfi.Size())
295 {
296 // Which boundary attributes need to be processed?
297 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
298 mesh->bdr_attributes.Max() : 0);
299 bdr_attr_marker = 0;
300 for (int k = 0; k < bnfi.Size(); k++)
301 {
302 if (bnfi_marker[k] == NULL)
303 {
304 bdr_attr_marker = 1;
305 break;
306 }
307 Array<int> &bdr_marker = *bnfi_marker[k];
308 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
309 "invalid boundary marker for boundary integrator #"
310 << k << ", counting from zero");
311 for (int i = 0; i < bdr_attr_marker.Size(); i++)
312 {
313 bdr_attr_marker[i] |= bdr_marker[i];
314 }
315 }
316
317 DofTransformation doftrans;
318 for (int i = 0; i < fes->GetNBE(); i++)
319 {
320 const int bdr_attr = mesh->GetBdrAttribute(i);
321 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
322
323 fe = fes->GetBE(i);
324
325 fes->GetBdrElementVDofs(i, vdofs, doftrans);
327 px.GetSubVector(vdofs, el_x);
328 doftrans.InvTransformPrimal(el_x);
329 for (int k = 0; k < bnfi.Size(); k++)
330 {
331 if (bnfi_marker[k] &&
332 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
333
334 bnfi[k]->AssembleElementVector(*fe, *T, el_x, el_y);
335 doftrans.TransformDual(el_y);
336 py.AddElementVector(vdofs, el_y);
337 }
338 }
339 }
340
341 if (fnfi.Size())
342 {
344 const FiniteElement *fe1, *fe2;
345 Array<int> vdofs2;
346
347 for (int i = 0; i < mesh->GetNumFaces(); i++)
348 {
349 tr = mesh->GetInteriorFaceTransformations(i);
350 if (tr != NULL)
351 {
352 fes->GetElementVDofs(tr->Elem1No, vdofs);
353 fes->GetElementVDofs(tr->Elem2No, vdofs2);
354 vdofs.Append (vdofs2);
355
356 px.GetSubVector(vdofs, el_x);
357
358 fe1 = fes->GetFE(tr->Elem1No);
359 fe2 = fes->GetFE(tr->Elem2No);
360
361 for (int k = 0; k < fnfi.Size(); k++)
362 {
363 fnfi[k]->AssembleFaceVector(*fe1, *fe2, *tr, el_x, el_y);
364 py.AddElementVector(vdofs, el_y);
365 }
366 }
367 }
368 }
369
370 if (bfnfi.Size())
371 {
373 const FiniteElement *fe1, *fe2;
374
375 // Which boundary attributes need to be processed?
376 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
377 mesh->bdr_attributes.Max() : 0);
378 bdr_attr_marker = 0;
379 for (int k = 0; k < bfnfi.Size(); k++)
380 {
381 if (bfnfi_marker[k] == NULL)
382 {
383 bdr_attr_marker = 1;
384 break;
385 }
386 Array<int> &bdr_marker = *bfnfi_marker[k];
387 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
388 "invalid boundary marker for boundary face integrator #"
389 << k << ", counting from zero");
390 for (int i = 0; i < bdr_attr_marker.Size(); i++)
391 {
392 bdr_attr_marker[i] |= bdr_marker[i];
393 }
394 }
395
396 for (int i = 0; i < fes -> GetNBE(); i++)
397 {
398 const int bdr_attr = mesh->GetBdrAttribute(i);
399 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
400
401 tr = mesh->GetBdrFaceTransformations (i);
402 if (tr != NULL)
403 {
404 fes->GetElementVDofs(tr->Elem1No, vdofs);
405 px.GetSubVector(vdofs, el_x);
406
407 fe1 = fes->GetFE(tr->Elem1No);
408 // The fe2 object is really a dummy and not used on the boundaries,
409 // but we can't dereference a NULL pointer, and we don't want to
410 // actually make a fake element.
411 fe2 = fe1;
412 for (int k = 0; k < bfnfi.Size(); k++)
413 {
414 if (bfnfi_marker[k] &&
415 (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
416
417 bfnfi[k]->AssembleFaceVector(*fe1, *fe2, *tr, el_x, el_y);
418 py.AddElementVector(vdofs, el_y);
419 }
420 }
421 }
422 }
423
424 if (Serial())
425 {
426 if (cP) { cP->MultTranspose(py, y); }
427
428 y.HostReadWrite();
430 for (int i = 0; i < ess_tdof_list.Size(); i++)
431 {
432 y(ess_tdof_list[i]) = 0.0;
433 }
434 // y(ess_tdof_list[i]) = x(ess_tdof_list[i]);
435 }
436 // In parallel, the result is in 'py' which is an alias for 'aux2'.
437}
438
439Operator &NonlinearForm::GetGradient(const Vector &x, bool finalize) const
440{
441 if (ext)
442 {
443 hGrad.Clear();
444 Operator &grad = ext->GetGradient(Prolongate(x));
445 Operator *Gop;
447 hGrad.Reset(Gop);
448 // In both serial and parallel, when using extension, we return the final
449 // global true-dof gradient with imposed b.c.
450 return *hGrad;
451 }
452
453 const int skip_zeros = 0;
454 Array<int> vdofs;
455 Vector el_x;
456 DenseMatrix elmat;
457 const FiniteElement *fe;
459 Mesh *mesh = fes->GetMesh();
460 const Vector &px = Prolongate(x);
461
462 if (Grad == NULL)
463 {
464 Grad = new SparseMatrix(fes->GetVSize());
465 }
466 else
467 {
468 *Grad = 0.0;
469 }
470
471 if (dnfi.Size())
472 {
473 // Which attributes need to be processed?
474 Array<int> attr_marker(mesh->attributes.Size() ?
475 mesh->attributes.Max() : 0);
476 attr_marker = 0;
477 for (int k = 0; k < dnfi.Size(); k++)
478 {
479 if (dnfi_marker[k] == NULL)
480 {
481 attr_marker = 1;
482 break;
483 }
484 Array<int> &marker = *dnfi_marker[k];
485 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
486 "invalid marker for domain integrator #"
487 << k << ", counting from zero");
488 for (int i = 0; i < attr_marker.Size(); i++)
489 {
490 attr_marker[i] |= marker[i];
491 }
492 }
493
494 DofTransformation doftrans;
495 for (int i = 0; i < fes->GetNE(); i++)
496 {
497 const int attr = mesh->GetAttribute(i);
498 if (attr_marker[attr-1] == 0) { continue; }
499
500 fe = fes->GetFE(i);
501 fes->GetElementVDofs(i, vdofs, doftrans);
503 px.GetSubVector(vdofs, el_x);
504 doftrans.InvTransformPrimal(el_x);
505 for (int k = 0; k < dnfi.Size(); k++)
506 {
507 if (dnfi_marker[k] &&
508 (*dnfi_marker[k])[attr-1] == 0) { continue; }
509
510 dnfi[k]->AssembleElementGrad(*fe, *T, el_x, elmat);
511 doftrans.TransformDual(elmat);
512 Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
513 // Grad->AddSubMatrix(vdofs, vdofs, elmat, 1);
514 }
515 }
516 }
517
518 if (bnfi.Size())
519 {
520 // Which boundary attributes need to be processed?
521 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
522 mesh->bdr_attributes.Max() : 0);
523 bdr_attr_marker = 0;
524 for (int k = 0; k < bnfi.Size(); k++)
525 {
526 if (bnfi_marker[k] == NULL)
527 {
528 bdr_attr_marker = 1;
529 break;
530 }
531 Array<int> &bdr_marker = *bnfi_marker[k];
532 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
533 "invalid boundary marker for boundary integrator #"
534 << k << ", counting from zero");
535 for (int i = 0; i < bdr_attr_marker.Size(); i++)
536 {
537 bdr_attr_marker[i] |= bdr_marker[i];
538 }
539 }
540
541 DofTransformation doftrans;
542 for (int i = 0; i < fes->GetNBE(); i++)
543 {
544 const int bdr_attr = mesh->GetBdrAttribute(i);
545 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
546
547 fe = fes->GetBE(i);
548 fes->GetBdrElementVDofs(i, vdofs, doftrans);
550 px.GetSubVector(vdofs, el_x);
551 doftrans.InvTransformPrimal(el_x);
552 for (int k = 0; k < bnfi.Size(); k++)
553 {
554 if (bnfi_marker[k] &&
555 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
556
557 bnfi[k]->AssembleElementGrad(*fe, *T, el_x, elmat);
558 doftrans.TransformDual(elmat);
559 Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
560 }
561 }
562 }
563
564 if (fnfi.Size())
565 {
567 const FiniteElement *fe1, *fe2;
568 Array<int> vdofs2;
569
570 for (int i = 0; i < mesh->GetNumFaces(); i++)
571 {
572 tr = mesh->GetInteriorFaceTransformations(i);
573 if (tr != NULL)
574 {
575 fes->GetElementVDofs(tr->Elem1No, vdofs);
576 fes->GetElementVDofs(tr->Elem2No, vdofs2);
577 vdofs.Append (vdofs2);
578
579 px.GetSubVector(vdofs, el_x);
580
581 fe1 = fes->GetFE(tr->Elem1No);
582 fe2 = fes->GetFE(tr->Elem2No);
583
584 for (int k = 0; k < fnfi.Size(); k++)
585 {
586 fnfi[k]->AssembleFaceGrad(*fe1, *fe2, *tr, el_x, elmat);
587 Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
588 }
589 }
590 }
591 }
592
593 if (bfnfi.Size())
594 {
596 const FiniteElement *fe1, *fe2;
597
598 // Which boundary attributes need to be processed?
599 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
600 mesh->bdr_attributes.Max() : 0);
601 bdr_attr_marker = 0;
602 for (int k = 0; k < bfnfi.Size(); k++)
603 {
604 if (bfnfi_marker[k] == NULL)
605 {
606 bdr_attr_marker = 1;
607 break;
608 }
609 Array<int> &bdr_marker = *bfnfi_marker[k];
610 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
611 "invalid boundary marker for boundary face integrator #"
612 << k << ", counting from zero");
613 for (int i = 0; i < bdr_attr_marker.Size(); i++)
614 {
615 bdr_attr_marker[i] |= bdr_marker[i];
616 }
617 }
618
619 for (int i = 0; i < fes -> GetNBE(); i++)
620 {
621 const int bdr_attr = mesh->GetBdrAttribute(i);
622 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
623
624 tr = mesh->GetBdrFaceTransformations (i);
625 if (tr != NULL)
626 {
627 fes->GetElementVDofs(tr->Elem1No, vdofs);
628 px.GetSubVector(vdofs, el_x);
629
630 fe1 = fes->GetFE(tr->Elem1No);
631 // The fe2 object is really a dummy and not used on the boundaries,
632 // but we can't dereference a NULL pointer, and we don't want to
633 // actually make a fake element.
634 fe2 = fe1;
635 for (int k = 0; k < bfnfi.Size(); k++)
636 {
637 if (bfnfi_marker[k] &&
638 (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
639
640 bfnfi[k]->AssembleFaceGrad(*fe1, *fe2, *tr, el_x, elmat);
641 Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
642 }
643 }
644 }
645 }
646
647 if (!finalize) { return *Grad; }
648
649 if (!Grad->Finalized())
650 {
651 Grad->Finalize(skip_zeros);
652 }
653
654 SparseMatrix *mGrad = Grad;
655 if (Serial())
656 {
657 if (cP)
658 {
659 delete cGrad;
660 cGrad = RAP(*cP, *Grad, *cP);
661 mGrad = cGrad;
662 }
663 for (int i = 0; i < ess_tdof_list.Size(); i++)
664 {
666 }
667 }
668
669 return *mGrad;
670}
671
673{
674 if (sequence == fes->GetSequence()) { return; }
675
677 delete cGrad; cGrad = NULL;
678 delete Grad; Grad = NULL;
679 hGrad.Clear();
680 ess_tdof_list.SetSize(0); // essential b.c. will need to be set again
682 // Do not modify aux1 and aux2, their size will be set before use.
684 cP = dynamic_cast<const SparseMatrix*>(P);
685
686 if (ext) { ext->Update(); }
687}
688
690{
691 if (ext) { ext->Assemble(); }
692}
693
695{
696 delete cGrad;
697 delete Grad;
698 if (!extern_bfs)
699 {
700 for (int i = 0; i < dnfi.Size(); i++) { delete dnfi[i]; }
701 for (int i = 0; i < bnfi.Size(); i++) { delete bnfi[i]; }
702 for (int i = 0; i < fnfi.Size(); i++) { delete fnfi[i]; }
703 for (int i = 0; i < bfnfi.Size(); i++) { delete bfnfi[i]; }
704 }
705 delete ext;
706}
707
708
710 fes(0), BlockGrad(NULL)
711{
712 height = 0;
713 width = 0;
714}
715
717{
718 delete BlockGrad;
719 BlockGrad = NULL;
720 for (int i=0; i<Grads.NumRows(); ++i)
721 {
722 for (int j=0; j<Grads.NumCols(); ++j)
723 {
724 delete Grads(i,j);
725 delete cGrads(i,j);
726 }
727 }
728 for (int i = 0; i < ess_tdofs.Size(); ++i)
729 {
730 delete ess_tdofs[i];
731 }
732
733 height = 0;
734 width = 0;
735 f.Copy(fes);
736 block_offsets.SetSize(f.Size() + 1);
737 block_trueOffsets.SetSize(f.Size() + 1);
738 block_offsets[0] = 0;
739 block_trueOffsets[0] = 0;
740
741 for (int i=0; i<fes.Size(); ++i)
742 {
743 block_offsets[i+1] = fes[i]->GetVSize();
744 block_trueOffsets[i+1] = fes[i]->GetTrueVSize();
745 }
746
749
750 height = block_trueOffsets[fes.Size()];
751 width = block_trueOffsets[fes.Size()];
752
753 Grads.SetSize(fes.Size(), fes.Size());
754 Grads = NULL;
755
756 cGrads.SetSize(fes.Size(), fes.Size());
757 cGrads = NULL;
758
759 P.SetSize(fes.Size());
760 cP.SetSize(fes.Size());
761 ess_tdofs.SetSize(fes.Size());
762 for (int s = 0; s < fes.Size(); ++s)
763 {
764 // Retrieve prolongation matrix for each FE space
765 P[s] = fes[s]->GetProlongationMatrix();
766 cP[s] = dynamic_cast<const SparseMatrix *>(P[s]);
767
768 // If the P Operator exists and its type is not SparseMatrix, this
769 // indicates the Operator is part of parallel run.
770 if (P[s] && !cP[s])
771 {
772 is_serial = false;
773 }
774
775 // If the P Operator exists and its type is SparseMatrix, this indicates
776 // the Operator is serial but needs prolongation on assembly.
777 if (cP[s])
778 {
779 needs_prolongation = true;
780 }
781
782 ess_tdofs[s] = new Array<int>;
783 }
784}
785
791
793 const Array<Array<int>*> &bdr_attr_is_ess, Array<Vector*> &rhs)
794{
795 for (int s = 0; s < fes.Size(); ++s)
796 {
797 fes[s]->GetEssentialTrueDofs(*bdr_attr_is_ess[s], *ess_tdofs[s]);
798
799 if (rhs[s])
800 {
801 rhs[s]->SetSubVector(*ess_tdofs[s], 0.0);
802 }
803 }
804}
805
807 const Array<Array<int>*> &ess_tdof_list, Array<Vector*> &rhs)
808{
809 for (int s = 0; s < fes.Size(); ++s)
810 {
811 *ess_tdofs[s] = *ess_tdof_list[s];
812 if (rhs[s])
813 {
814 rhs[s]->SetSubVector(*ess_tdofs[s], 0.0);
815 }
816 }
817}
818
820{
821 Array<Array<int> *> vdofs(fes.Size());
822 Array<Vector *> el_x(fes.Size());
823 Array<const Vector *> el_x_const(fes.Size());
826 Mesh *mesh = fes[0]->GetMesh();
827 real_t energy = 0.0;
828
829 for (int i=0; i<fes.Size(); ++i)
830 {
831 el_x_const[i] = el_x[i] = new Vector();
832 vdofs[i] = new Array<int>;
833 }
834
835 if (dnfi.Size())
836 {
837 // Which attributes need to be processed?
838 Array<int> attr_marker(mesh->attributes.Size() ?
839 mesh->attributes.Max() : 0);
840 attr_marker = 0;
841 for (int k = 0; k < dnfi.Size(); k++)
842 {
843 if (dnfi_marker[k] == NULL)
844 {
845 attr_marker = 1;
846 break;
847 }
848 Array<int> &marker = *dnfi_marker[k];
849 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
850 "invalid marker for domain integrator #"
851 << k << ", counting from zero");
852 for (int i = 0; i < attr_marker.Size(); i++)
853 {
854 attr_marker[i] |= marker[i];
855 }
856 }
857
858 DofTransformation doftrans;
859 for (int i = 0; i < fes[0]->GetNE(); ++i)
860 {
861 const int attr = mesh->GetAttribute(i);
862 if (attr_marker[attr-1] == 0) { continue; }
863
864 T = fes[0]->GetElementTransformation(i);
865 for (int s=0; s<fes.Size(); ++s)
866 {
867 fe[s] = fes[s]->GetFE(i);
868 fes[s]->GetElementVDofs(i, *vdofs[s], doftrans);
869 bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
870 doftrans.InvTransformPrimal(*el_x[s]);
871 }
872
873 for (int k = 0; k < dnfi.Size(); ++k)
874 {
875 if (dnfi_marker[k] &&
876 (*dnfi_marker[k])[attr-1] == 0) { continue; }
877
878 energy += dnfi[k]->GetElementEnergy(fe, *T, el_x_const);
879 }
880 }
881 }
882
883 if (bnfi.Size())
884 {
885 // Which boundary attributes need to be processed?
886 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
887 mesh->bdr_attributes.Max() : 0);
888 bdr_attr_marker = 0;
889 for (int k = 0; k < bnfi.Size(); k++)
890 {
891 if (bnfi_marker[k] == NULL)
892 {
893 bdr_attr_marker = 1;
894 break;
895 }
896 Array<int> &bdr_marker = *bnfi_marker[k];
897 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
898 "invalid boundary marker for boundary integrator #"
899 << k << ", counting from zero");
900 for (int i = 0; i < bdr_attr_marker.Size(); i++)
901 {
902 bdr_attr_marker[i] |= bdr_marker[i];
903 }
904 }
905
906 DofTransformation doftrans;
907 for (int i = 0; i < mesh->GetNBE(); i++)
908 {
909 const int bdr_attr = mesh->GetBdrAttribute(i);
910 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
911
912 T = fes[0]->GetBdrElementTransformation(i);
913 for (int s = 0; s < fes.Size(); ++s)
914 {
915 fe[s] = fes[s]->GetBE(i);
916 fes[s]->GetBdrElementVDofs(i, *(vdofs[s]), doftrans);
917 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
918 doftrans.InvTransformPrimal(*el_x[s]);
919 }
920
921 for (int k = 0; k < bnfi.Size(); k++)
922 {
923 if (bnfi_marker[k] &&
924 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
925
926 energy += bnfi[k]->GetElementEnergy(fe, *T, el_x_const);
927 }
928 }
929 }
930
931 // free the allocated memory
932 for (int i = 0; i < fes.Size(); ++i)
933 {
934 delete el_x[i];
935 delete vdofs[i];
936 }
937
938 if (fnfi.Size())
939 {
940 MFEM_ABORT("TODO: add energy contribution from interior face terms");
941 }
942
943 if (bfnfi.Size())
944 {
945 MFEM_ABORT("TODO: add energy contribution from boundary face terms");
946 }
947
948 return energy;
949}
950
952{
953 xs.Update(const_cast<Vector&>(x), block_offsets);
954 return GetEnergyBlocked(xs);
955}
956
958 BlockVector &by) const
959{
960 Array<Array<int> *>vdofs(fes.Size());
961 Array<Array<int> *>vdofs2(fes.Size());
962 Array<Vector *> el_x(fes.Size());
963 Array<const Vector *> el_x_const(fes.Size());
964 Array<Vector *> el_y(fes.Size());
968 std::vector<DofTransformation> doftrans(fes.Size());
969 Mesh *mesh = fes[0]->GetMesh();
970
971 by.UseDevice(true);
972 by = 0.0;
973 by.SyncToBlocks();
974 for (int s=0; s<fes.Size(); ++s)
975 {
976 el_x_const[s] = el_x[s] = new Vector();
977 el_y[s] = new Vector();
978 vdofs[s] = new Array<int>;
979 vdofs2[s] = new Array<int>;
980 }
981
982 if (dnfi.Size())
983 {
984 // Which attributes need to be processed?
985 Array<int> attr_marker(mesh->attributes.Size() ?
986 mesh->attributes.Max() : 0);
987 attr_marker = 0;
988 for (int k = 0; k < dnfi.Size(); k++)
989 {
990 if (dnfi_marker[k] == NULL)
991 {
992 attr_marker = 1;
993 break;
994 }
995 Array<int> &marker = *dnfi_marker[k];
996 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
997 "invalid marker for domain integrator #"
998 << k << ", counting from zero");
999 for (int i = 0; i < attr_marker.Size(); i++)
1000 {
1001 attr_marker[i] |= marker[i];
1002 }
1003 }
1004
1005 for (int i = 0; i < fes[0]->GetNE(); ++i)
1006 {
1007 const int attr = mesh->GetAttribute(i);
1008 if (attr_marker[attr-1] == 0) { continue; }
1009
1010 T = fes[0]->GetElementTransformation(i);
1011 for (int s = 0; s < fes.Size(); ++s)
1012 {
1013 fes[s]->GetElementVDofs(i, *(vdofs[s]), doftrans[s]);
1014 fe[s] = fes[s]->GetFE(i);
1015 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
1016 doftrans[s].InvTransformPrimal(*el_x[s]);
1017 }
1018
1019 for (int k = 0; k < dnfi.Size(); ++k)
1020 {
1021 if (dnfi_marker[k] &&
1022 (*dnfi_marker[k])[attr-1] == 0) { continue; }
1023
1024 dnfi[k]->AssembleElementVector(fe, *T,
1025 el_x_const, el_y);
1026
1027 for (int s=0; s<fes.Size(); ++s)
1028 {
1029 if (el_y[s]->Size() == 0) { continue; }
1030 doftrans[s].TransformDual(*el_y[s]);
1031 by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
1032 }
1033 }
1034 }
1035 }
1036
1037 if (bnfi.Size())
1038 {
1039 // Which boundary attributes need to be processed?
1040 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1041 mesh->bdr_attributes.Max() : 0);
1042 bdr_attr_marker = 0;
1043 for (int k = 0; k < bnfi.Size(); k++)
1044 {
1045 if (bnfi_marker[k] == NULL)
1046 {
1047 bdr_attr_marker = 1;
1048 break;
1049 }
1050 Array<int> &bdr_marker = *bnfi_marker[k];
1051 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1052 "invalid boundary marker for boundary integrator #"
1053 << k << ", counting from zero");
1054 for (int i = 0; i < bdr_attr_marker.Size(); i++)
1055 {
1056 bdr_attr_marker[i] |= bdr_marker[i];
1057 }
1058 }
1059
1060 for (int i = 0; i < mesh->GetNBE(); i++)
1061 {
1062 const int bdr_attr = mesh->GetBdrAttribute(i);
1063 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1064
1065 T = fes[0]->GetBdrElementTransformation(i);
1066 for (int s = 0; s < fes.Size(); ++s)
1067 {
1068 fes[s]->GetBdrElementVDofs(i, *(vdofs[s]), doftrans[s]);
1069 fe[s] = fes[s]->GetBE(i);
1070 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
1071 doftrans[s].InvTransformPrimal(*el_x[s]);
1072 }
1073
1074 for (int k = 0; k < bnfi.Size(); k++)
1075 {
1076 if (bnfi_marker[k] &&
1077 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
1078
1079 bnfi[k]->AssembleElementVector(fe, *T, el_x_const, el_y);
1080
1081 for (int s=0; s<fes.Size(); ++s)
1082 {
1083 if (el_y[s]->Size() == 0) { continue; }
1084 doftrans[s].TransformDual(*el_y[s]);
1085 by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
1086 }
1087 }
1088 }
1089 }
1090
1091
1092 if (fnfi.Size())
1093 {
1095
1096 for (int i = 0; i < mesh->GetNumFaces(); ++i)
1097 {
1098 tr = mesh->GetInteriorFaceTransformations(i);
1099 if (tr != NULL)
1100 {
1101 for (int s=0; s<fes.Size(); ++s)
1102 {
1103 fe[s] = fes[s]->GetFE(tr->Elem1No);
1104 fe2[s] = fes[s]->GetFE(tr->Elem2No);
1105
1106 fes[s]->GetElementVDofs(tr->Elem1No, *(vdofs[s]));
1107 fes[s]->GetElementVDofs(tr->Elem2No, *(vdofs2[s]));
1108
1109 vdofs[s]->Append(*(vdofs2[s]));
1110
1111 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
1112 }
1113
1114 for (int k = 0; k < fnfi.Size(); ++k)
1115 {
1116
1117 fnfi[k]->AssembleFaceVector(fe, fe2, *tr, el_x_const, el_y);
1118
1119 for (int s=0; s<fes.Size(); ++s)
1120 {
1121 if (el_y[s]->Size() == 0) { continue; }
1122 by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
1123 }
1124 }
1125 }
1126 }
1127 }
1128
1129 if (bfnfi.Size())
1130 {
1132
1133 // Which boundary attributes need to be processed?
1134 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1135 mesh->bdr_attributes.Max() : 0);
1136 bdr_attr_marker = 0;
1137 for (int k = 0; k < bfnfi.Size(); ++k)
1138 {
1139 if (bfnfi_marker[k] == NULL)
1140 {
1141 bdr_attr_marker = 1;
1142 break;
1143 }
1144 Array<int> &bdr_marker = *bfnfi_marker[k];
1145 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1146 "invalid boundary marker for boundary face integrator #"
1147 << k << ", counting from zero");
1148 for (int i = 0; i < bdr_attr_marker.Size(); ++i)
1149 {
1150 bdr_attr_marker[i] |= bdr_marker[i];
1151 }
1152 }
1153
1154 for (int i = 0; i < mesh->GetNBE(); ++i)
1155 {
1156 const int bdr_attr = mesh->GetBdrAttribute(i);
1157 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1158
1159 tr = mesh->GetBdrFaceTransformations(i);
1160 if (tr != NULL)
1161 {
1162 for (int s=0; s<fes.Size(); ++s)
1163 {
1164 fe[s] = fes[s]->GetFE(tr->Elem1No);
1165 fe2[s] = fes[s]->GetFE(tr->Elem1No);
1166
1167 fes[s]->GetElementVDofs(tr->Elem1No, *(vdofs[s]));
1168 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
1169 }
1170
1171 for (int k = 0; k < bfnfi.Size(); ++k)
1172 {
1173 if (bfnfi_marker[k] &&
1174 (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
1175
1176 bfnfi[k]->AssembleFaceVector(fe, fe2, *tr, el_x_const, el_y);
1177
1178 for (int s=0; s<fes.Size(); ++s)
1179 {
1180 if (el_y[s]->Size() == 0) { continue; }
1181 by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
1182 }
1183 }
1184 }
1185 }
1186 }
1187
1188 for (int s=0; s<fes.Size(); ++s)
1189 {
1190 delete vdofs2[s];
1191 delete vdofs[s];
1192 delete el_y[s];
1193 delete el_x[s];
1194 }
1195
1196 by.SyncFromBlocks();
1197}
1198
1200{
1201 MFEM_VERIFY(bx.Size() == Width(), "invalid input BlockVector size");
1202
1204 {
1206 for (int s = 0; s < fes.Size(); s++)
1207 {
1208 if (P[s])
1209 {
1210 P[s]->Mult(bx.GetBlock(s), aux1.GetBlock(s));
1211 }
1212 else
1213 {
1214 aux1.GetBlock(s) = bx.GetBlock(s);
1215 }
1216 }
1217 return aux1;
1218 }
1219 return bx;
1220}
1221
1223{
1224 BlockVector bx(const_cast<Vector&>(x), block_trueOffsets);
1226
1227 const BlockVector &pbx = Prolongate(bx);
1229 {
1231 }
1232 BlockVector &pby = needs_prolongation ? aux2 : by;
1233
1234 xs.Update(const_cast<BlockVector&>(pbx), block_offsets);
1235 ys.Update(pby, block_offsets);
1236 MultBlocked(xs, ys);
1237
1238 for (int s = 0; s < fes.Size(); s++)
1239 {
1240 if (cP[s])
1241 {
1242 cP[s]->MultTranspose(pby.GetBlock(s), by.GetBlock(s));
1243 }
1244 else if (needs_prolongation)
1245 {
1246 by.GetBlock(s) = pby.GetBlock(s);
1247 }
1248 by.GetBlock(s).SetSubVector(*ess_tdofs[s], 0.0);
1249 }
1250}
1251
1253 bool finalize) const
1254{
1255 const int skip_zeros = 0;
1256 Array<Array<int> *> vdofs(fes.Size());
1257 Array<Array<int> *> vdofs2(fes.Size());
1258 Array<Vector *> el_x(fes.Size());
1259 Array<const Vector *> el_x_const(fes.Size());
1260 Array2D<DenseMatrix *> elmats(fes.Size(), fes.Size());
1264 std::vector<DofTransformation> doftrans(fes.Size());
1265 Mesh *mesh = fes[0]->GetMesh();
1266
1267 for (int i=0; i<fes.Size(); ++i)
1268 {
1269 el_x_const[i] = el_x[i] = new Vector();
1270 vdofs[i] = new Array<int>;
1271 vdofs2[i] = new Array<int>;
1272 for (int j=0; j<fes.Size(); ++j)
1273 {
1274 elmats(i,j) = new DenseMatrix();
1275 }
1276 }
1277
1278 for (int i=0; i<fes.Size(); ++i)
1279 {
1280 for (int j=0; j<fes.Size(); ++j)
1281 {
1282 if (Grads(i,j) != NULL)
1283 {
1284 *Grads(i,j) = 0.0;
1285 }
1286 else
1287 {
1288 Grads(i,j) = new SparseMatrix(fes[i]->GetVSize(),
1289 fes[j]->GetVSize());
1290 }
1291 }
1292 }
1293
1294 if (dnfi.Size())
1295 {
1296 // Which attributes need to be processed?
1297 Array<int> attr_marker(mesh->attributes.Size() ?
1298 mesh->attributes.Max() : 0);
1299 attr_marker = 0;
1300 for (int k = 0; k < dnfi.Size(); k++)
1301 {
1302 if (dnfi_marker[k] == NULL)
1303 {
1304 attr_marker = 1;
1305 break;
1306 }
1307 Array<int> &marker = *dnfi_marker[k];
1308 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
1309 "invalid marker for domain integrator #"
1310 << k << ", counting from zero");
1311 for (int i = 0; i < attr_marker.Size(); i++)
1312 {
1313 attr_marker[i] |= marker[i];
1314 }
1315 }
1316
1317 for (int i = 0; i < fes[0]->GetNE(); ++i)
1318 {
1319 const int attr = mesh->GetAttribute(i);
1320 if (attr_marker[attr-1] == 0) { continue; }
1321
1322 T = fes[0]->GetElementTransformation(i);
1323 for (int s = 0; s < fes.Size(); ++s)
1324 {
1325 fe[s] = fes[s]->GetFE(i);
1326 fes[s]->GetElementVDofs(i, *vdofs[s], doftrans[s]);
1327 bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
1328 doftrans[s].InvTransformPrimal(*el_x[s]);
1329 }
1330
1331 for (int k = 0; k < dnfi.Size(); ++k)
1332 {
1333 if (dnfi_marker[k] &&
1334 (*dnfi_marker[k])[attr-1] == 0) { continue; }
1335
1336 dnfi[k]->AssembleElementGrad(fe, *T, el_x_const, elmats);
1337
1338 for (int j=0; j<fes.Size(); ++j)
1339 {
1340 for (int l=0; l<fes.Size(); ++l)
1341 {
1342 if (elmats(j,l)->Height() == 0) { continue; }
1343 TransformDual(doftrans[j], doftrans[l], *elmats(j,l));
1344 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1345 *elmats(j,l), skip_zeros);
1346 }
1347 }
1348 }
1349 }
1350 }
1351
1352 if (bnfi.Size())
1353 {
1354 // Which boundary attributes need to be processed?
1355 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1356 mesh->bdr_attributes.Max() : 0);
1357 bdr_attr_marker = 0;
1358 for (int k = 0; k < bnfi.Size(); k++)
1359 {
1360 if (bnfi_marker[k] == NULL)
1361 {
1362 bdr_attr_marker = 1;
1363 break;
1364 }
1365 Array<int> &bdr_marker = *bnfi_marker[k];
1366 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1367 "invalid boundary marker for boundary integrator #"
1368 << k << ", counting from zero");
1369 for (int i = 0; i < bdr_attr_marker.Size(); i++)
1370 {
1371 bdr_attr_marker[i] |= bdr_marker[i];
1372 }
1373 }
1374
1375 for (int i = 0; i < mesh->GetNBE(); i++)
1376 {
1377 const int bdr_attr = mesh->GetBdrAttribute(i);
1378 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1379
1380 T = fes[0]->GetBdrElementTransformation(i);
1381 for (int s = 0; s < fes.Size(); ++s)
1382 {
1383 fe[s] = fes[s]->GetBE(i);
1384 fes[s]->GetBdrElementVDofs(i, *(vdofs[s]), doftrans[s]);
1385 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
1386 doftrans[s].InvTransformPrimal(*el_x[s]);
1387 }
1388
1389 for (int k = 0; k < bnfi.Size(); k++)
1390 {
1391 if (bnfi_marker[k] &&
1392 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
1393
1394 bnfi[k]->AssembleElementGrad(fe, *T, el_x_const, elmats);
1395
1396 for (int j=0; j<fes.Size(); ++j)
1397 {
1398 for (int l=0; l<fes.Size(); ++l)
1399 {
1400 if (elmats(j,l)->Height() == 0) { continue; }
1401 TransformDual(doftrans[j], doftrans[l], *elmats(j,l));
1402 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1403 *elmats(j,l), skip_zeros);
1404 }
1405 }
1406 }
1407 }
1408 }
1409
1410 if (fnfi.Size())
1411 {
1413
1414 for (int i = 0; i < mesh->GetNumFaces(); ++i)
1415 {
1416 tr = mesh->GetInteriorFaceTransformations(i);
1417
1418 for (int s=0; s < fes.Size(); ++s)
1419 {
1420 fe[s] = fes[s]->GetFE(tr->Elem1No);
1421 fe2[s] = fes[s]->GetFE(tr->Elem2No);
1422
1423 fes[s]->GetElementVDofs(tr->Elem1No, *vdofs[s]);
1424 fes[s]->GetElementVDofs(tr->Elem2No, *vdofs2[s]);
1425 vdofs[s]->Append(*(vdofs2[s]));
1426
1427 bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
1428 }
1429
1430 for (int k = 0; k < fnfi.Size(); ++k)
1431 {
1432 fnfi[k]->AssembleFaceGrad(fe, fe2, *tr, el_x_const, elmats);
1433 for (int j=0; j<fes.Size(); ++j)
1434 {
1435 for (int l=0; l<fes.Size(); ++l)
1436 {
1437 if (elmats(j,l)->Height() == 0) { continue; }
1438 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1439 *elmats(j,l), skip_zeros);
1440 }
1441 }
1442 }
1443 }
1444 }
1445
1446 if (bfnfi.Size())
1447 {
1449
1450 // Which boundary attributes need to be processed?
1451 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1452 mesh->bdr_attributes.Max() : 0);
1453 bdr_attr_marker = 0;
1454 for (int k = 0; k < bfnfi.Size(); ++k)
1455 {
1456 if (bfnfi_marker[k] == NULL)
1457 {
1458 bdr_attr_marker = 1;
1459 break;
1460 }
1461 Array<int> &bdr_marker = *bfnfi_marker[k];
1462 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1463 "invalid boundary marker for boundary face integrator #"
1464 << k << ", counting from zero");
1465 for (int i = 0; i < bdr_attr_marker.Size(); ++i)
1466 {
1467 bdr_attr_marker[i] |= bdr_marker[i];
1468 }
1469 }
1470
1471 for (int i = 0; i < mesh->GetNBE(); ++i)
1472 {
1473 const int bdr_attr = mesh->GetBdrAttribute(i);
1474 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1475
1476 tr = mesh->GetBdrFaceTransformations(i);
1477 if (tr != NULL)
1478 {
1479 for (int s = 0; s < fes.Size(); ++s)
1480 {
1481 fe[s] = fes[s]->GetFE(tr->Elem1No);
1482 fe2[s] = fe[s];
1483
1484 fes[s]->GetElementVDofs(tr->Elem1No, *vdofs[s]);
1485 bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
1486 }
1487
1488 for (int k = 0; k < bfnfi.Size(); ++k)
1489 {
1490 if (bfnfi_marker[k] &&
1491 (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
1492 bfnfi[k]->AssembleFaceGrad(fe, fe2, *tr, el_x_const, elmats);
1493 for (int l=0; l<fes.Size(); ++l)
1494 {
1495 for (int j=0; j<fes.Size(); ++j)
1496 {
1497 if (elmats(j,l)->Height() == 0) { continue; }
1498 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1499 *elmats(j,l), skip_zeros);
1500 }
1501 }
1502 }
1503 }
1504 }
1505 }
1506
1507 if (finalize && !Grads(0,0)->Finalized())
1508 {
1509 for (int i=0; i<fes.Size(); ++i)
1510 {
1511 for (int j=0; j<fes.Size(); ++j)
1512 {
1513 Grads(i,j)->Finalize(skip_zeros);
1514 }
1515 }
1516 }
1517
1518 for (int i=0; i<fes.Size(); ++i)
1519 {
1520 for (int j=0; j<fes.Size(); ++j)
1521 {
1522 delete elmats(i,j);
1523 }
1524 delete vdofs2[i];
1525 delete vdofs[i];
1526 delete el_x[i];
1527 }
1528}
1529
1531{
1532 BlockVector bx(const_cast<Vector&>(x), block_trueOffsets);
1533 const BlockVector &pbx = Prolongate(bx);
1534
1536
1537 Array2D<SparseMatrix *> mGrads(fes.Size(), fes.Size());
1538 mGrads = Grads;
1540 {
1541 for (int s1 = 0; s1 < fes.Size(); ++s1)
1542 {
1543 for (int s2 = 0; s2 < fes.Size(); ++s2)
1544 {
1545 delete cGrads(s1, s2);
1546 if (cP[s1] && cP[s2])
1547 {
1548 cGrads(s1, s2) = RAP(*cP[s1], *Grads(s1, s2), *cP[s2]);
1549 }
1550 else if (cP[s1])
1551 {
1552 cGrads(s1, s2) = TransposeMult(*cP[s1], *Grads(s1, s2));
1553 }
1554 else if (cP[s2])
1555 {
1556 cGrads(s1, s2) = mfem::Mult(*Grads(s1, s2), *cP[s2]);
1557 }
1558 else
1559 {
1560 cGrads(s1, s2) = NULL;
1561 continue;
1562 }
1563 mGrads(s1, s2) = cGrads(s1, s2);
1564 }
1565 }
1566 }
1567
1568 for (int s = 0; s < fes.Size(); ++s)
1569 {
1570 for (int i = 0; i < ess_tdofs[s]->Size(); ++i)
1571 {
1572 for (int j = 0; j < fes.Size(); ++j)
1573 {
1574 if (s == j)
1575 {
1576 mGrads(s, s)->EliminateRowCol((*ess_tdofs[s])[i],
1578 }
1579 else
1580 {
1581 mGrads(s, j)->EliminateRow((*ess_tdofs[s])[i]);
1582 mGrads(j, s)->EliminateCol((*ess_tdofs[s])[i]);
1583 }
1584 }
1585 }
1586 }
1587
1588 delete BlockGrad;
1590 for (int i = 0; i < fes.Size(); ++i)
1591 {
1592 for (int j = 0; j < fes.Size(); ++j)
1593 {
1594 BlockGrad->SetBlock(i, j, mGrads(i, j));
1595 }
1596 }
1597 return *BlockGrad;
1598}
1599
1601{
1602 delete BlockGrad;
1603 for (int i=0; i<fes.Size(); ++i)
1604 {
1605 for (int j=0; j<fes.Size(); ++j)
1606 {
1607 delete Grads(i,j);
1608 delete cGrads(i,j);
1609 }
1610 delete ess_tdofs[i];
1611 }
1612
1613 for (int i = 0; i < dnfi.Size(); ++i)
1614 {
1615 delete dnfi[i];
1616 }
1617
1618 for (int i = 0; i < bnfi.Size(); ++i)
1619 {
1620 delete bnfi[i];
1621 }
1622
1623 for (int i = 0; i < fnfi.Size(); ++i)
1624 {
1625 delete fnfi[i];
1626 }
1627
1628 for (int i = 0; i < bfnfi.Size(); ++i)
1629 {
1630 delete bfnfi[i];
1631 }
1632
1633}
1634
1635}
Dynamic 2D array using row-major layout.
Definition array.hpp:430
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:69
const T * HostRead() const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), false).
Definition array.hpp:385
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
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:912
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition array.hpp:1042
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:381
Array< Array< int > * > dnfi_marker
Array< Array< int > * > bnfi_marker
Array< BlockNonlinearFormIntegrator * > dnfi
Set of Domain Integrators to be assembled (added).
void Mult(const Vector &x, Vector &y) const override
Operator & GetGradient(const Vector &x) const override
Array< FiniteElementSpace * > fes
FE spaces on which the form lives.
real_t GetEnergyBlocked(const BlockVector &bx) const
Specialized version of GetEnergy() for BlockVectors.
void ComputeGradientBlocked(const BlockVector &bx, bool finalize=true) const
Specialized version of GetGradient() for BlockVector.
void SetSpaces(Array< FiniteElementSpace * > &f)
(Re)initialize the BlockNonlinearForm.
Array2D< SparseMatrix * > Grads
bool is_serial
Indicator if the Operator is part of a parallel run.
void MultBlocked(const BlockVector &bx, BlockVector &by) const
virtual ~BlockNonlinearForm()
Destructor.
Array< BlockNonlinearFormIntegrator * > fnfi
Set of interior face Integrators to be assembled (added).
const BlockVector & Prolongate(const BlockVector &bx) const
Array< Array< int > * > bfnfi_marker
virtual void SetEssentialBC(const Array< Array< int > * > &bdr_attr_is_ess, Array< Vector * > &rhs)
Set essential boundary conditions to each finite element space using boundary attribute markers.
bool needs_prolongation
Indicator if the Operator needs prolongation on assembly.
Array2D< SparseMatrix * > cGrads
Array< BlockNonlinearFormIntegrator * > bnfi
Set of Boundary Integrators to be assembled (added).
Array< const SparseMatrix * > cP
Array of results of dynamic-casting P to SparseMatrix pointer.
Array< BlockNonlinearFormIntegrator * > bfnfi
Set of Boundary Face Integrators to be assembled (added).
virtual void SetEssentialTrueDofs(const Array< Array< int > * > &ess_tdof_list, Array< Vector * > &rhs)
Set essential boundary conditions to each finite element space using essential true dof lists.
virtual real_t GetEnergy(const Vector &x) const
Array< const Operator * > P
Array of pointers to the prolongation matrix of fes, may be NULL.
Array< Array< int > * > ess_tdofs
BlockNonlinearForm()
Construct an empty BlockNonlinearForm. Initialize with SetSpaces().
A class to handle Block systems in a matrix-free implementation.
void SetBlock(int iRow, int iCol, Operator *op, real_t c=1.0)
Add a block op in the block-entry (iblock, jblock).
A class to handle Vectors in a block fashion.
void Update(real_t *data, const Array< int > &bOffsets)
Update method.
void SyncFromBlocks() const
Synchronize the memory location flags (i.e. the memory validity flags) of the big/monolithic block-ve...
void SyncToBlocks() const
Synchronize the memory location flags (i.e. the memory validity flags) of the big/monolithic block-ve...
Vector & GetBlock(int i)
Get the i-th vector in the block.
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void TransformDual(real_t *v) const
Definition doftrans.cpp:77
void InvTransformPrimal(real_t *v) const
Definition doftrans.cpp:47
A specialized ElementTransformation class representing a face and its two neighboring elements.
Definition eltrans.hpp:750
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th boundary fac...
Definition fespace.cpp:3870
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:829
ElementTransformation * GetElementTransformation(int i) const
Definition fespace.hpp:905
static void ListToMarker(const Array< int > &list, int marker_size, Array< int > &marker, int mark_val=-1)
Convert an array of indices (list) to a Boolean marker array where all indices in the list are marked...
Definition fespace.cpp:779
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition fespace.cpp:628
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition fespace.hpp:878
virtual const Operator * GetProlongationMatrix() const
Definition fespace.hpp:696
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
Definition fespace.cpp:304
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition fespace.cpp:3824
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i-th boundary element.
Definition fespace.hpp:914
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:869
static void MarkerToList(const Array< int > &marker, Array< int > &list)
Convert a Boolean marker array to a list containing all marked indices.
Definition fespace.cpp:760
void ConvertToConformingVDofs(const Array< int > &dofs, Array< int > &cdofs)
For a partially conforming FE space, convert a marker array (nonzero entries are true) on the partial...
Definition fespace.cpp:792
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:644
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:826
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for i'th boundary element. The returned indices are offsets int...
Definition fespace.cpp:319
Abstract class for all finite elements.
Definition fe_base.hpp:247
void BooleanMultTranspose(int alpha, const int *x, int beta, int *y)
The "Boolean" analog of y = alpha * A^T * x + beta * y, where elements in the sparsity pattern of the...
Definition hypre.hpp:817
Data and methods for unassembled nonlinear forms.
Mesh data type.
Definition mesh.hpp:65
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:304
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition mesh.cpp:6846
int GetAttribute(int i) const
Return the attribute of element i.
Definition mesh.hpp:1484
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition mesh.hpp:1490
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Builds the transformation defining the given boundary face.
Definition mesh.cpp:1223
int GetNBE() const
Returns number of boundary elements.
Definition mesh.hpp:1380
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
See GetFaceElementTransformations().
Definition mesh.cpp:1203
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition mesh.hpp:302
virtual Operator & GetGradient(const Vector &x) const =0
Return the gradient as an L-to-L Operator. The input x must be an L-vector (i.e. GridFunction-size ve...
virtual void Update()=0
Called by NonlinearForm::Update() to reflect changes in the FE space.
virtual void Assemble()=0
Assemble at the AssemblyLevel of the subclass.
virtual real_t GetGridFunctionEnergy(const Vector &x) const =0
Compute the local (to the MPI rank) energy of the L-vector state x.
Array< Array< int > * > dnfi_marker
FiniteElementSpace * fes
FE space on which the form lives.
long sequence
Counter for updates propagated from the FiniteElementSpace.
bool extern_bfs
Indicates the NonlinearFormIntegrators stored in dnfi, bnfi, fnfi, and bfnfi are not owned by this No...
AssemblyLevel assembly
The assembly level.
SparseMatrix * Grad
const Vector & Prolongate(const Vector &x) const
Vector aux1
Auxiliary Vectors.
const SparseMatrix * cP
The result of dynamic-casting P to SparseMatrix pointer.
NonlinearFormExtension * ext
Extension for supporting different AssemblyLevels.
Operator & GetGradient(const Vector &x) const override
Compute the gradient Operator of the NonlinearForm corresponding to the state x.
Array< NonlinearFormIntegrator * > dnfi
Set of Domain Integrators to be assembled (added).
OperatorHandle hGrad
Gradient Operator when not assembled as a matrix.
const Operator * P
Pointer to the prolongation matrix of fes, may be NULL.
Array< NonlinearFormIntegrator * > bnfi
Set of Boundary Integrators to be assembled (added).
virtual void Update()
Update the NonlinearForm to propagate updates of the associated FE space.
Array< int > ess_tdof_list
A list of all essential true dofs.
Array< Array< int > * > bfnfi_marker
virtual ~NonlinearForm()
Destroy the NonlinearForm including the owned NonlinearFormIntegrators and gradient Operator.
SparseMatrix * cGrad
Array< Array< int > * > bnfi_marker
void Mult(const Vector &x, Vector &y) const override
Evaluate the action of the NonlinearForm.
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACY.
void SetEssentialBC(const Array< int > &bdr_attr_is_ess, Vector *rhs=NULL)
Specify essential boundary conditions.
Array< NonlinearFormIntegrator * > fnfi
Set of interior face Integrators to be assembled (added).
Array< NonlinearFormIntegrator * > bfnfi
Set of boundary face Integrators to be assembled (added).
real_t GetGridFunctionEnergy(const Vector &x) const
Compute the energy corresponding to the state x.
virtual void Setup()
Setup the NonlinearForm: based on the current AssemblyLevel and the current mesh, optionally,...
void SetEssentialVDofs(const Array< int > &ess_vdofs_list)
Specify essential boundary conditions.
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged.
Definition handle.hpp:124
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition handle.hpp:145
Abstract operator.
Definition operator.hpp:25
int width
Dimension of the input / number of columns in the matrix.
Definition operator.hpp:28
void FormSystemOperator(const Array< int > &ess_tdof_list, Operator *&A)
Return in A a parallel (on truedofs) version of this square operator.
Definition operator.cpp:227
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
int height
Dimension of the output / number of rows in the matrix.
Definition operator.hpp:27
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
@ DIAG_ONE
Set the diagonal value to one.
Definition operator.hpp:50
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition operator.hpp:72
Data and methods for partially-assembled nonlinear forms.
Abstract parallel finite element space.
Definition pfespace.hpp:31
int GetTrueVSize() const override
Return the number of local vector true dofs.
Definition pfespace.hpp:353
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition pfespace.hpp:391
Data type sparse matrix.
Definition sparsemat.hpp:51
void EliminateRowCol(int rc, const real_t sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Eliminate row rc and column rc and modify the rhs using sol.
void MultTranspose(const Vector &x, Vector &y) const override
Multiply a vector with the transposed matrix. y = At * x.
bool Finalized() const
Returns whether or not CSR format has been finalized.
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
void Finalize(int skip_zeros=1) override
Finalize the matrix initialization, switching the storage format from LIL to CSR.
Vector data type.
Definition vector.hpp:82
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:536
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
Definition vector.cpp:702
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
int Size() const
Returns the size of the vector.
Definition vector.hpp:234
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition vector.hpp:145
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:584
virtual real_t * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:540
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition vector.cpp:676
void mfem_error(const char *msg)
Definition error.cpp:154
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
AssemblyLevel
Enumeration defining the assembly level for bilinear and nonlinear form classes derived from Operator...
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
float real_t
Definition config.hpp:46
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
SparseMatrix * TransposeMult(const SparseMatrix &A, const SparseMatrix &B)
C = A^T B.
void forall(int N, lambda &&body)
Definition forall.hpp:839