MFEM v4.8.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 DofTransformation *doftrans;
100 Mesh *mesh = fes->GetMesh();
101 real_t energy = 0.0;
102
103 if (dnfi.Size())
104 {
105 // Which attributes need to be processed?
106 Array<int> attr_marker(mesh->attributes.Size() ?
107 mesh->attributes.Max() : 0);
108 attr_marker = 0;
109 for (int k = 0; k < dnfi.Size(); k++)
110 {
111 if (dnfi_marker[k] == NULL)
112 {
113 attr_marker = 1;
114 break;
115 }
116 Array<int> &marker = *dnfi_marker[k];
117 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
118 "invalid marker for domain integrator #"
119 << k << ", counting from zero");
120 for (int i = 0; i < attr_marker.Size(); i++)
121 {
122 attr_marker[i] |= marker[i];
123 }
124 }
125
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 doftrans = fes->GetElementVDofs(i, vdofs);
134 x.GetSubVector(vdofs, el_x);
135 if (doftrans) {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 for (int i = 0; i < fes->GetNBE(); i++)
170 {
171 const int bdr_attr = mesh->GetBdrAttribute(i);
172 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
173
174 fe = fes->GetBE(i);
175 doftrans = fes->GetBdrElementVDofs(i, vdofs);
177 x.GetSubVector(vdofs, el_x);
178 if (doftrans) {doftrans->InvTransformPrimal(el_x); }
179 for (int k = 0; k < bnfi.Size(); k++)
180 {
181 if (bnfi_marker[k] &&
182 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
183
184 energy += bnfi[k]->GetElementEnergy(*fe, *T, el_x);
185 }
186 }
187
188 }
189
190 if (fnfi.Size())
191 {
192 MFEM_ABORT("TODO: add energy contribution from interior face terms");
193 }
194
195 if (bfnfi.Size())
196 {
197 MFEM_ABORT("TODO: add energy contribution from boundary face terms");
198 }
199
200 return energy;
201}
202
204{
205 MFEM_VERIFY(x.Size() == Width(), "invalid input Vector size");
206 if (P)
207 {
208 aux1.SetSize(P->Height());
209 P->Mult(x, aux1);
210 return aux1;
211 }
212 return x;
213}
214
215void NonlinearForm::Mult(const Vector &x, Vector &y) const
216{
217 const Vector &px = Prolongate(x);
218 if (P) { aux2.SetSize(P->Height()); }
219
220 // If we are in parallel, ParNonLinearForm::Mult uses the aux2 vector. In
221 // serial, place the result directly in y (when there is no P).
222 Vector &py = P ? aux2 : y;
223
224 if (ext)
225 {
226 ext->Mult(px, py);
227 if (Serial())
228 {
229 if (cP) { cP->MultTranspose(py, y); }
230 const int N = ess_tdof_list.Size();
231 const auto tdof = ess_tdof_list.Read();
232 auto Y = y.ReadWrite();
233 mfem::forall(N, [=] MFEM_HOST_DEVICE (int i) { Y[tdof[i]] = 0.0; });
234 }
235 // In parallel, the result is in 'py' which is an alias for 'aux2'.
236 return;
237 }
238
239 Array<int> vdofs;
240 Vector el_x, el_y;
241 const FiniteElement *fe;
243 DofTransformation *doftrans;
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 for (int i = 0; i < fes->GetNE(); i++)
272 {
273 const int attr = mesh->GetAttribute(i);
274 if (attr_marker[attr-1] == 0) { continue; }
275
276 fe = fes->GetFE(i);
277 doftrans = fes->GetElementVDofs(i, vdofs);
279 px.GetSubVector(vdofs, el_x);
280 if (doftrans) {doftrans->InvTransformPrimal(el_x); }
281 for (int k = 0; k < dnfi.Size(); k++)
282 {
283 if (dnfi_marker[k] &&
284 (*dnfi_marker[k])[attr-1] == 0) { continue; }
285
286 dnfi[k]->AssembleElementVector(*fe, *T, el_x, el_y);
287 if (doftrans) {doftrans->TransformDual(el_y); }
288 py.AddElementVector(vdofs, el_y);
289 }
290 }
291 }
292
293 if (bnfi.Size())
294 {
295 // Which boundary attributes need to be processed?
296 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
297 mesh->bdr_attributes.Max() : 0);
298 bdr_attr_marker = 0;
299 for (int k = 0; k < bnfi.Size(); k++)
300 {
301 if (bnfi_marker[k] == NULL)
302 {
303 bdr_attr_marker = 1;
304 break;
305 }
306 Array<int> &bdr_marker = *bnfi_marker[k];
307 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
308 "invalid boundary marker for boundary integrator #"
309 << k << ", counting from zero");
310 for (int i = 0; i < bdr_attr_marker.Size(); i++)
311 {
312 bdr_attr_marker[i] |= bdr_marker[i];
313 }
314 }
315
316 for (int i = 0; i < fes->GetNBE(); i++)
317 {
318 const int bdr_attr = mesh->GetBdrAttribute(i);
319 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
320
321 fe = fes->GetBE(i);
322 doftrans = fes->GetBdrElementVDofs(i, vdofs);
324 px.GetSubVector(vdofs, el_x);
325 if (doftrans) {doftrans->InvTransformPrimal(el_x); }
326 for (int k = 0; k < bnfi.Size(); k++)
327 {
328 if (bnfi_marker[k] &&
329 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
330
331 bnfi[k]->AssembleElementVector(*fe, *T, el_x, el_y);
332 if (doftrans) {doftrans->TransformDual(el_y); }
333 py.AddElementVector(vdofs, el_y);
334 }
335 }
336 }
337
338 if (fnfi.Size())
339 {
341 const FiniteElement *fe1, *fe2;
342 Array<int> vdofs2;
343
344 for (int i = 0; i < mesh->GetNumFaces(); i++)
345 {
346 tr = mesh->GetInteriorFaceTransformations(i);
347 if (tr != NULL)
348 {
349 fes->GetElementVDofs(tr->Elem1No, vdofs);
350 fes->GetElementVDofs(tr->Elem2No, vdofs2);
351 vdofs.Append (vdofs2);
352
353 px.GetSubVector(vdofs, el_x);
354
355 fe1 = fes->GetFE(tr->Elem1No);
356 fe2 = fes->GetFE(tr->Elem2No);
357
358 for (int k = 0; k < fnfi.Size(); k++)
359 {
360 fnfi[k]->AssembleFaceVector(*fe1, *fe2, *tr, el_x, el_y);
361 py.AddElementVector(vdofs, el_y);
362 }
363 }
364 }
365 }
366
367 if (bfnfi.Size())
368 {
370 const FiniteElement *fe1, *fe2;
371
372 // Which boundary attributes need to be processed?
373 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
374 mesh->bdr_attributes.Max() : 0);
375 bdr_attr_marker = 0;
376 for (int k = 0; k < bfnfi.Size(); k++)
377 {
378 if (bfnfi_marker[k] == NULL)
379 {
380 bdr_attr_marker = 1;
381 break;
382 }
383 Array<int> &bdr_marker = *bfnfi_marker[k];
384 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
385 "invalid boundary marker for boundary face integrator #"
386 << k << ", counting from zero");
387 for (int i = 0; i < bdr_attr_marker.Size(); i++)
388 {
389 bdr_attr_marker[i] |= bdr_marker[i];
390 }
391 }
392
393 for (int i = 0; i < fes -> GetNBE(); i++)
394 {
395 const int bdr_attr = mesh->GetBdrAttribute(i);
396 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
397
398 tr = mesh->GetBdrFaceTransformations (i);
399 if (tr != NULL)
400 {
401 fes->GetElementVDofs(tr->Elem1No, vdofs);
402 px.GetSubVector(vdofs, el_x);
403
404 fe1 = fes->GetFE(tr->Elem1No);
405 // The fe2 object is really a dummy and not used on the boundaries,
406 // but we can't dereference a NULL pointer, and we don't want to
407 // actually make a fake element.
408 fe2 = fe1;
409 for (int k = 0; k < bfnfi.Size(); k++)
410 {
411 if (bfnfi_marker[k] &&
412 (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
413
414 bfnfi[k]->AssembleFaceVector(*fe1, *fe2, *tr, el_x, el_y);
415 py.AddElementVector(vdofs, el_y);
416 }
417 }
418 }
419 }
420
421 if (Serial())
422 {
423 if (cP) { cP->MultTranspose(py, y); }
424
425 y.HostReadWrite();
427 for (int i = 0; i < ess_tdof_list.Size(); i++)
428 {
429 y(ess_tdof_list[i]) = 0.0;
430 }
431 // y(ess_tdof_list[i]) = x(ess_tdof_list[i]);
432 }
433 // In parallel, the result is in 'py' which is an alias for 'aux2'.
434}
435
437{
438 if (ext)
439 {
440 hGrad.Clear();
441 Operator &grad = ext->GetGradient(Prolongate(x));
442 Operator *Gop;
444 hGrad.Reset(Gop);
445 // In both serial and parallel, when using extension, we return the final
446 // global true-dof gradient with imposed b.c.
447 return *hGrad;
448 }
449
450 const int skip_zeros = 0;
451 Array<int> vdofs;
452 Vector el_x;
453 DenseMatrix elmat;
454 const FiniteElement *fe;
456 DofTransformation *doftrans;
457 Mesh *mesh = fes->GetMesh();
458 const Vector &px = Prolongate(x);
459
460 if (Grad == NULL)
461 {
462 Grad = new SparseMatrix(fes->GetVSize());
463 }
464 else
465 {
466 *Grad = 0.0;
467 }
468
469 if (dnfi.Size())
470 {
471 // Which attributes need to be processed?
472 Array<int> attr_marker(mesh->attributes.Size() ?
473 mesh->attributes.Max() : 0);
474 attr_marker = 0;
475 for (int k = 0; k < dnfi.Size(); k++)
476 {
477 if (dnfi_marker[k] == NULL)
478 {
479 attr_marker = 1;
480 break;
481 }
482 Array<int> &marker = *dnfi_marker[k];
483 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
484 "invalid marker for domain integrator #"
485 << k << ", counting from zero");
486 for (int i = 0; i < attr_marker.Size(); i++)
487 {
488 attr_marker[i] |= marker[i];
489 }
490 }
491
492 for (int i = 0; i < fes->GetNE(); i++)
493 {
494 const int attr = mesh->GetAttribute(i);
495 if (attr_marker[attr-1] == 0) { continue; }
496
497 fe = fes->GetFE(i);
498 doftrans = fes->GetElementVDofs(i, vdofs);
500 px.GetSubVector(vdofs, el_x);
501 if (doftrans) {doftrans->InvTransformPrimal(el_x); }
502 for (int k = 0; k < dnfi.Size(); k++)
503 {
504 if (dnfi_marker[k] &&
505 (*dnfi_marker[k])[attr-1] == 0) { continue; }
506
507 dnfi[k]->AssembleElementGrad(*fe, *T, el_x, elmat);
508 if (doftrans) { doftrans->TransformDual(elmat); }
509 Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
510 // Grad->AddSubMatrix(vdofs, vdofs, elmat, 1);
511 }
512 }
513 }
514
515 if (bnfi.Size())
516 {
517 // Which boundary attributes need to be processed?
518 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
519 mesh->bdr_attributes.Max() : 0);
520 bdr_attr_marker = 0;
521 for (int k = 0; k < bnfi.Size(); k++)
522 {
523 if (bnfi_marker[k] == NULL)
524 {
525 bdr_attr_marker = 1;
526 break;
527 }
528 Array<int> &bdr_marker = *bnfi_marker[k];
529 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
530 "invalid boundary marker for boundary integrator #"
531 << k << ", counting from zero");
532 for (int i = 0; i < bdr_attr_marker.Size(); i++)
533 {
534 bdr_attr_marker[i] |= bdr_marker[i];
535 }
536 }
537
538 for (int i = 0; i < fes->GetNBE(); i++)
539 {
540 const int bdr_attr = mesh->GetBdrAttribute(i);
541 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
542
543 fe = fes->GetBE(i);
544 doftrans = fes->GetBdrElementVDofs(i, vdofs);
546 px.GetSubVector(vdofs, el_x);
547 if (doftrans) {doftrans->InvTransformPrimal(el_x); }
548 for (int k = 0; k < bnfi.Size(); k++)
549 {
550 if (bnfi_marker[k] &&
551 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
552
553 bnfi[k]->AssembleElementGrad(*fe, *T, el_x, elmat);
554 if (doftrans) { doftrans->TransformDual(elmat); }
555 Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
556 }
557 }
558 }
559
560 if (fnfi.Size())
561 {
563 const FiniteElement *fe1, *fe2;
564 Array<int> vdofs2;
565
566 for (int i = 0; i < mesh->GetNumFaces(); i++)
567 {
568 tr = mesh->GetInteriorFaceTransformations(i);
569 if (tr != NULL)
570 {
571 fes->GetElementVDofs(tr->Elem1No, vdofs);
572 fes->GetElementVDofs(tr->Elem2No, vdofs2);
573 vdofs.Append (vdofs2);
574
575 px.GetSubVector(vdofs, el_x);
576
577 fe1 = fes->GetFE(tr->Elem1No);
578 fe2 = fes->GetFE(tr->Elem2No);
579
580 for (int k = 0; k < fnfi.Size(); k++)
581 {
582 fnfi[k]->AssembleFaceGrad(*fe1, *fe2, *tr, el_x, elmat);
583 Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
584 }
585 }
586 }
587 }
588
589 if (bfnfi.Size())
590 {
592 const FiniteElement *fe1, *fe2;
593
594 // Which boundary attributes need to be processed?
595 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
596 mesh->bdr_attributes.Max() : 0);
597 bdr_attr_marker = 0;
598 for (int k = 0; k < bfnfi.Size(); k++)
599 {
600 if (bfnfi_marker[k] == NULL)
601 {
602 bdr_attr_marker = 1;
603 break;
604 }
605 Array<int> &bdr_marker = *bfnfi_marker[k];
606 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
607 "invalid boundary marker for boundary face integrator #"
608 << k << ", counting from zero");
609 for (int i = 0; i < bdr_attr_marker.Size(); i++)
610 {
611 bdr_attr_marker[i] |= bdr_marker[i];
612 }
613 }
614
615 for (int i = 0; i < fes -> GetNBE(); i++)
616 {
617 const int bdr_attr = mesh->GetBdrAttribute(i);
618 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
619
620 tr = mesh->GetBdrFaceTransformations (i);
621 if (tr != NULL)
622 {
623 fes->GetElementVDofs(tr->Elem1No, vdofs);
624 px.GetSubVector(vdofs, el_x);
625
626 fe1 = fes->GetFE(tr->Elem1No);
627 // The fe2 object is really a dummy and not used on the boundaries,
628 // but we can't dereference a NULL pointer, and we don't want to
629 // actually make a fake element.
630 fe2 = fe1;
631 for (int k = 0; k < bfnfi.Size(); k++)
632 {
633 if (bfnfi_marker[k] &&
634 (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
635
636 bfnfi[k]->AssembleFaceGrad(*fe1, *fe2, *tr, el_x, elmat);
637 Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
638 }
639 }
640 }
641 }
642
643 if (!Grad->Finalized())
644 {
645 Grad->Finalize(skip_zeros);
646 }
647
648 SparseMatrix *mGrad = Grad;
649 if (Serial())
650 {
651 if (cP)
652 {
653 delete cGrad;
654 cGrad = RAP(*cP, *Grad, *cP);
655 mGrad = cGrad;
656 }
657 for (int i = 0; i < ess_tdof_list.Size(); i++)
658 {
660 }
661 }
662
663 return *mGrad;
664}
665
667{
668 if (sequence == fes->GetSequence()) { return; }
669
671 delete cGrad; cGrad = NULL;
672 delete Grad; Grad = NULL;
673 hGrad.Clear();
674 ess_tdof_list.SetSize(0); // essential b.c. will need to be set again
676 // Do not modify aux1 and aux2, their size will be set before use.
678 cP = dynamic_cast<const SparseMatrix*>(P);
679
680 if (ext) { ext->Update(); }
681}
682
684{
685 if (ext) { ext->Assemble(); }
686}
687
689{
690 delete cGrad;
691 delete Grad;
692 if (!extern_bfs)
693 {
694 for (int i = 0; i < dnfi.Size(); i++) { delete dnfi[i]; }
695 for (int i = 0; i < bnfi.Size(); i++) { delete bnfi[i]; }
696 for (int i = 0; i < fnfi.Size(); i++) { delete fnfi[i]; }
697 for (int i = 0; i < bfnfi.Size(); i++) { delete bfnfi[i]; }
698 }
699 delete ext;
700}
701
702
704 fes(0), BlockGrad(NULL)
705{
706 height = 0;
707 width = 0;
708}
709
711{
712 delete BlockGrad;
713 BlockGrad = NULL;
714 for (int i=0; i<Grads.NumRows(); ++i)
715 {
716 for (int j=0; j<Grads.NumCols(); ++j)
717 {
718 delete Grads(i,j);
719 delete cGrads(i,j);
720 }
721 }
722 for (int i = 0; i < ess_tdofs.Size(); ++i)
723 {
724 delete ess_tdofs[i];
725 }
726
727 height = 0;
728 width = 0;
729 f.Copy(fes);
730 block_offsets.SetSize(f.Size() + 1);
731 block_trueOffsets.SetSize(f.Size() + 1);
732 block_offsets[0] = 0;
733 block_trueOffsets[0] = 0;
734
735 for (int i=0; i<fes.Size(); ++i)
736 {
737 block_offsets[i+1] = fes[i]->GetVSize();
738 block_trueOffsets[i+1] = fes[i]->GetTrueVSize();
739 }
740
743
744 height = block_trueOffsets[fes.Size()];
745 width = block_trueOffsets[fes.Size()];
746
747 Grads.SetSize(fes.Size(), fes.Size());
748 Grads = NULL;
749
750 cGrads.SetSize(fes.Size(), fes.Size());
751 cGrads = NULL;
752
753 P.SetSize(fes.Size());
754 cP.SetSize(fes.Size());
755 ess_tdofs.SetSize(fes.Size());
756 for (int s = 0; s < fes.Size(); ++s)
757 {
758 // Retrieve prolongation matrix for each FE space
759 P[s] = fes[s]->GetProlongationMatrix();
760 cP[s] = dynamic_cast<const SparseMatrix *>(P[s]);
761
762 // If the P Operator exists and its type is not SparseMatrix, this
763 // indicates the Operator is part of parallel run.
764 if (P[s] && !cP[s])
765 {
766 is_serial = false;
767 }
768
769 // If the P Operator exists and its type is SparseMatrix, this indicates
770 // the Operator is serial but needs prolongation on assembly.
771 if (cP[s])
772 {
773 needs_prolongation = true;
774 }
775
776 ess_tdofs[s] = new Array<int>;
777 }
778}
779
785
787 const Array<Array<int> *> &bdr_attr_is_ess, Array<Vector *> &rhs)
788{
789 for (int s = 0; s < fes.Size(); ++s)
790 {
791 ess_tdofs[s]->SetSize(ess_tdofs.Size());
792
793 fes[s]->GetEssentialTrueDofs(*bdr_attr_is_ess[s], *ess_tdofs[s]);
794
795 if (rhs[s])
796 {
797 rhs[s]->SetSubVector(*ess_tdofs[s], 0.0);
798 }
799 }
800}
801
803{
804 Array<Array<int> *> vdofs(fes.Size());
805 Array<Vector *> el_x(fes.Size());
806 Array<const Vector *> el_x_const(fes.Size());
809 DofTransformation *doftrans;
810 Mesh *mesh = fes[0]->GetMesh();
811 real_t energy = 0.0;
812
813 for (int i=0; i<fes.Size(); ++i)
814 {
815 el_x_const[i] = el_x[i] = new Vector();
816 vdofs[i] = new Array<int>;
817 }
818
819 if (dnfi.Size())
820 {
821 // Which attributes need to be processed?
822 Array<int> attr_marker(mesh->attributes.Size() ?
823 mesh->attributes.Max() : 0);
824 attr_marker = 0;
825 for (int k = 0; k < dnfi.Size(); k++)
826 {
827 if (dnfi_marker[k] == NULL)
828 {
829 attr_marker = 1;
830 break;
831 }
832 Array<int> &marker = *dnfi_marker[k];
833 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
834 "invalid marker for domain integrator #"
835 << k << ", counting from zero");
836 for (int i = 0; i < attr_marker.Size(); i++)
837 {
838 attr_marker[i] |= marker[i];
839 }
840 }
841
842 for (int i = 0; i < fes[0]->GetNE(); ++i)
843 {
844 const int attr = mesh->GetAttribute(i);
845 if (attr_marker[attr-1] == 0) { continue; }
846
847 T = fes[0]->GetElementTransformation(i);
848 for (int s=0; s<fes.Size(); ++s)
849 {
850 fe[s] = fes[s]->GetFE(i);
851 doftrans = fes[s]->GetElementVDofs(i, *vdofs[s]);
852 bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
853 if (doftrans) {doftrans->InvTransformPrimal(*el_x[s]); }
854 }
855
856 for (int k = 0; k < dnfi.Size(); ++k)
857 {
858 if (dnfi_marker[k] &&
859 (*dnfi_marker[k])[attr-1] == 0) { continue; }
860
861 energy += dnfi[k]->GetElementEnergy(fe, *T, el_x_const);
862 }
863 }
864 }
865
866 if (bnfi.Size())
867 {
868 // Which boundary attributes need to be processed?
869 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
870 mesh->bdr_attributes.Max() : 0);
871 bdr_attr_marker = 0;
872 for (int k = 0; k < bnfi.Size(); k++)
873 {
874 if (bnfi_marker[k] == NULL)
875 {
876 bdr_attr_marker = 1;
877 break;
878 }
879 Array<int> &bdr_marker = *bnfi_marker[k];
880 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
881 "invalid boundary marker for boundary integrator #"
882 << k << ", counting from zero");
883 for (int i = 0; i < bdr_attr_marker.Size(); i++)
884 {
885 bdr_attr_marker[i] |= bdr_marker[i];
886 }
887 }
888
889 for (int i = 0; i < mesh->GetNBE(); i++)
890 {
891 const int bdr_attr = mesh->GetBdrAttribute(i);
892 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
893
894 T = fes[0]->GetBdrElementTransformation(i);
895 for (int s = 0; s < fes.Size(); ++s)
896 {
897 fe[s] = fes[s]->GetBE(i);
898 doftrans = fes[s]->GetBdrElementVDofs(i, *(vdofs[s]));
899 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
900 if (doftrans) {doftrans->InvTransformPrimal(*el_x[s]); }
901 }
902
903 for (int k = 0; k < bnfi.Size(); k++)
904 {
905 if (bnfi_marker[k] &&
906 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
907
908 energy += bnfi[k]->GetElementEnergy(fe, *T, el_x_const);
909 }
910 }
911 }
912
913 // free the allocated memory
914 for (int i = 0; i < fes.Size(); ++i)
915 {
916 delete el_x[i];
917 delete vdofs[i];
918 }
919
920 if (fnfi.Size())
921 {
922 MFEM_ABORT("TODO: add energy contribution from interior face terms");
923 }
924
925 if (bfnfi.Size())
926 {
927 MFEM_ABORT("TODO: add energy contribution from boundary face terms");
928 }
929
930 return energy;
931}
932
934{
935 xs.Update(const_cast<Vector&>(x), block_offsets);
936 return GetEnergyBlocked(xs);
937}
938
940 BlockVector &by) const
941{
942 Array<Array<int> *>vdofs(fes.Size());
943 Array<Array<int> *>vdofs2(fes.Size());
944 Array<Vector *> el_x(fes.Size());
945 Array<const Vector *> el_x_const(fes.Size());
946 Array<Vector *> el_y(fes.Size());
950 Array<DofTransformation *> doftrans(fes.Size()); doftrans = nullptr;
951 Mesh *mesh = fes[0]->GetMesh();
952
953 by.UseDevice(true);
954 by = 0.0;
955 by.SyncToBlocks();
956 for (int s=0; s<fes.Size(); ++s)
957 {
958 el_x_const[s] = el_x[s] = new Vector();
959 el_y[s] = new Vector();
960 vdofs[s] = new Array<int>;
961 vdofs2[s] = new Array<int>;
962 }
963
964 if (dnfi.Size())
965 {
966 // Which attributes need to be processed?
967 Array<int> attr_marker(mesh->attributes.Size() ?
968 mesh->attributes.Max() : 0);
969 attr_marker = 0;
970 for (int k = 0; k < dnfi.Size(); k++)
971 {
972 if (dnfi_marker[k] == NULL)
973 {
974 attr_marker = 1;
975 break;
976 }
977 Array<int> &marker = *dnfi_marker[k];
978 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
979 "invalid marker for domain integrator #"
980 << k << ", counting from zero");
981 for (int i = 0; i < attr_marker.Size(); i++)
982 {
983 attr_marker[i] |= marker[i];
984 }
985 }
986
987 for (int i = 0; i < fes[0]->GetNE(); ++i)
988 {
989 const int attr = mesh->GetAttribute(i);
990 if (attr_marker[attr-1] == 0) { continue; }
991
992 T = fes[0]->GetElementTransformation(i);
993 for (int s = 0; s < fes.Size(); ++s)
994 {
995 doftrans[s] = fes[s]->GetElementVDofs(i, *(vdofs[s]));
996 fe[s] = fes[s]->GetFE(i);
997 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
998 if (doftrans[s]) {doftrans[s]->InvTransformPrimal(*el_x[s]); }
999 }
1000
1001 for (int k = 0; k < dnfi.Size(); ++k)
1002 {
1003 if (dnfi_marker[k] &&
1004 (*dnfi_marker[k])[attr-1] == 0) { continue; }
1005
1006 dnfi[k]->AssembleElementVector(fe, *T,
1007 el_x_const, el_y);
1008
1009 for (int s=0; s<fes.Size(); ++s)
1010 {
1011 if (el_y[s]->Size() == 0) { continue; }
1012 if (doftrans[s]) {doftrans[s]->TransformDual(*el_y[s]); }
1013 by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
1014 }
1015 }
1016 }
1017 }
1018
1019 if (bnfi.Size())
1020 {
1021 // Which boundary attributes need to be processed?
1022 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1023 mesh->bdr_attributes.Max() : 0);
1024 bdr_attr_marker = 0;
1025 for (int k = 0; k < bnfi.Size(); k++)
1026 {
1027 if (bnfi_marker[k] == NULL)
1028 {
1029 bdr_attr_marker = 1;
1030 break;
1031 }
1032 Array<int> &bdr_marker = *bnfi_marker[k];
1033 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1034 "invalid boundary marker for boundary integrator #"
1035 << k << ", counting from zero");
1036 for (int i = 0; i < bdr_attr_marker.Size(); i++)
1037 {
1038 bdr_attr_marker[i] |= bdr_marker[i];
1039 }
1040 }
1041
1042 for (int i = 0; i < mesh->GetNBE(); i++)
1043 {
1044 const int bdr_attr = mesh->GetBdrAttribute(i);
1045 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1046
1047 T = fes[0]->GetBdrElementTransformation(i);
1048 for (int s = 0; s < fes.Size(); ++s)
1049 {
1050 doftrans[s] = fes[s]->GetBdrElementVDofs(i, *(vdofs[s]));
1051 fe[s] = fes[s]->GetBE(i);
1052 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
1053 if (doftrans[s]) {doftrans[s]->InvTransformPrimal(*el_x[s]); }
1054 }
1055
1056 for (int k = 0; k < bnfi.Size(); k++)
1057 {
1058 if (bnfi_marker[k] &&
1059 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
1060
1061 bnfi[k]->AssembleElementVector(fe, *T, el_x_const, el_y);
1062
1063 for (int s=0; s<fes.Size(); ++s)
1064 {
1065 if (el_y[s]->Size() == 0) { continue; }
1066 if (doftrans[s]) {doftrans[s]->TransformDual(*el_y[s]); }
1067 by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
1068 }
1069 }
1070 }
1071 }
1072
1073
1074 if (fnfi.Size())
1075 {
1077
1078 for (int i = 0; i < mesh->GetNumFaces(); ++i)
1079 {
1080 tr = mesh->GetInteriorFaceTransformations(i);
1081 if (tr != NULL)
1082 {
1083 for (int s=0; s<fes.Size(); ++s)
1084 {
1085 fe[s] = fes[s]->GetFE(tr->Elem1No);
1086 fe2[s] = fes[s]->GetFE(tr->Elem2No);
1087
1088 fes[s]->GetElementVDofs(tr->Elem1No, *(vdofs[s]));
1089 fes[s]->GetElementVDofs(tr->Elem2No, *(vdofs2[s]));
1090
1091 vdofs[s]->Append(*(vdofs2[s]));
1092
1093 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
1094 }
1095
1096 for (int k = 0; k < fnfi.Size(); ++k)
1097 {
1098
1099 fnfi[k]->AssembleFaceVector(fe, fe2, *tr, el_x_const, el_y);
1100
1101 for (int s=0; s<fes.Size(); ++s)
1102 {
1103 if (el_y[s]->Size() == 0) { continue; }
1104 by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
1105 }
1106 }
1107 }
1108 }
1109 }
1110
1111 if (bfnfi.Size())
1112 {
1114
1115 // Which boundary attributes need to be processed?
1116 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1117 mesh->bdr_attributes.Max() : 0);
1118 bdr_attr_marker = 0;
1119 for (int k = 0; k < bfnfi.Size(); ++k)
1120 {
1121 if (bfnfi_marker[k] == NULL)
1122 {
1123 bdr_attr_marker = 1;
1124 break;
1125 }
1126 Array<int> &bdr_marker = *bfnfi_marker[k];
1127 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1128 "invalid boundary marker for boundary face integrator #"
1129 << k << ", counting from zero");
1130 for (int i = 0; i < bdr_attr_marker.Size(); ++i)
1131 {
1132 bdr_attr_marker[i] |= bdr_marker[i];
1133 }
1134 }
1135
1136 for (int i = 0; i < mesh->GetNBE(); ++i)
1137 {
1138 const int bdr_attr = mesh->GetBdrAttribute(i);
1139 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1140
1141 tr = mesh->GetBdrFaceTransformations(i);
1142 if (tr != NULL)
1143 {
1144 for (int s=0; s<fes.Size(); ++s)
1145 {
1146 fe[s] = fes[s]->GetFE(tr->Elem1No);
1147 fe2[s] = fes[s]->GetFE(tr->Elem1No);
1148
1149 fes[s]->GetElementVDofs(tr->Elem1No, *(vdofs[s]));
1150 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
1151 }
1152
1153 for (int k = 0; k < bfnfi.Size(); ++k)
1154 {
1155 if (bfnfi_marker[k] &&
1156 (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
1157
1158 bfnfi[k]->AssembleFaceVector(fe, fe2, *tr, el_x_const, el_y);
1159
1160 for (int s=0; s<fes.Size(); ++s)
1161 {
1162 if (el_y[s]->Size() == 0) { continue; }
1163 by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
1164 }
1165 }
1166 }
1167 }
1168 }
1169
1170 for (int s=0; s<fes.Size(); ++s)
1171 {
1172 delete vdofs2[s];
1173 delete vdofs[s];
1174 delete el_y[s];
1175 delete el_x[s];
1176 }
1177
1178 by.SyncFromBlocks();
1179}
1180
1182{
1183 MFEM_VERIFY(bx.Size() == Width(), "invalid input BlockVector size");
1184
1186 {
1188 for (int s = 0; s < fes.Size(); s++)
1189 {
1190 P[s]->Mult(bx.GetBlock(s), aux1.GetBlock(s));
1191 }
1192 return aux1;
1193 }
1194 return bx;
1195}
1196
1198{
1199 BlockVector bx(const_cast<Vector&>(x), block_trueOffsets);
1201
1202 const BlockVector &pbx = Prolongate(bx);
1204 {
1206 }
1207 BlockVector &pby = needs_prolongation ? aux2 : by;
1208
1209 xs.Update(const_cast<BlockVector&>(pbx), block_offsets);
1210 ys.Update(pby, block_offsets);
1211 MultBlocked(xs, ys);
1212
1213 for (int s = 0; s < fes.Size(); s++)
1214 {
1215 if (cP[s])
1216 {
1217 cP[s]->MultTranspose(pby.GetBlock(s), by.GetBlock(s));
1218 }
1219 by.GetBlock(s).SetSubVector(*ess_tdofs[s], 0.0);
1220 }
1221}
1222
1224{
1225 const int skip_zeros = 0;
1226 Array<Array<int> *> vdofs(fes.Size());
1227 Array<Array<int> *> vdofs2(fes.Size());
1228 Array<Vector *> el_x(fes.Size());
1229 Array<const Vector *> el_x_const(fes.Size());
1230 Array2D<DenseMatrix *> elmats(fes.Size(), fes.Size());
1234 Array<DofTransformation *> doftrans(fes.Size()); doftrans = nullptr;
1235 Mesh *mesh = fes[0]->GetMesh();
1236
1237 for (int i=0; i<fes.Size(); ++i)
1238 {
1239 el_x_const[i] = el_x[i] = new Vector();
1240 vdofs[i] = new Array<int>;
1241 vdofs2[i] = new Array<int>;
1242 for (int j=0; j<fes.Size(); ++j)
1243 {
1244 elmats(i,j) = new DenseMatrix();
1245 }
1246 }
1247
1248 for (int i=0; i<fes.Size(); ++i)
1249 {
1250 for (int j=0; j<fes.Size(); ++j)
1251 {
1252 if (Grads(i,j) != NULL)
1253 {
1254 *Grads(i,j) = 0.0;
1255 }
1256 else
1257 {
1258 Grads(i,j) = new SparseMatrix(fes[i]->GetVSize(),
1259 fes[j]->GetVSize());
1260 }
1261 }
1262 }
1263
1264 if (dnfi.Size())
1265 {
1266 // Which attributes need to be processed?
1267 Array<int> attr_marker(mesh->attributes.Size() ?
1268 mesh->attributes.Max() : 0);
1269 attr_marker = 0;
1270 for (int k = 0; k < dnfi.Size(); k++)
1271 {
1272 if (dnfi_marker[k] == NULL)
1273 {
1274 attr_marker = 1;
1275 break;
1276 }
1277 Array<int> &marker = *dnfi_marker[k];
1278 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
1279 "invalid marker for domain integrator #"
1280 << k << ", counting from zero");
1281 for (int i = 0; i < attr_marker.Size(); i++)
1282 {
1283 attr_marker[i] |= marker[i];
1284 }
1285 }
1286
1287 for (int i = 0; i < fes[0]->GetNE(); ++i)
1288 {
1289 const int attr = mesh->GetAttribute(i);
1290 if (attr_marker[attr-1] == 0) { continue; }
1291
1292 T = fes[0]->GetElementTransformation(i);
1293 for (int s = 0; s < fes.Size(); ++s)
1294 {
1295 fe[s] = fes[s]->GetFE(i);
1296 doftrans[s] = fes[s]->GetElementVDofs(i, *vdofs[s]);
1297 bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
1298 if (doftrans[s]) {doftrans[s]->InvTransformPrimal(*el_x[s]); }
1299 }
1300
1301 for (int k = 0; k < dnfi.Size(); ++k)
1302 {
1303 if (dnfi_marker[k] &&
1304 (*dnfi_marker[k])[attr-1] == 0) { continue; }
1305
1306 dnfi[k]->AssembleElementGrad(fe, *T, el_x_const, elmats);
1307
1308 for (int j=0; j<fes.Size(); ++j)
1309 {
1310 for (int l=0; l<fes.Size(); ++l)
1311 {
1312 if (elmats(j,l)->Height() == 0) { continue; }
1313 if (doftrans[j] || doftrans[l])
1314 {
1315 TransformDual(doftrans[j], doftrans[l], *elmats(j,l));
1316 }
1317 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1318 *elmats(j,l), skip_zeros);
1319 }
1320 }
1321 }
1322 }
1323 }
1324
1325 if (bnfi.Size())
1326 {
1327 // Which boundary attributes need to be processed?
1328 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1329 mesh->bdr_attributes.Max() : 0);
1330 bdr_attr_marker = 0;
1331 for (int k = 0; k < bnfi.Size(); k++)
1332 {
1333 if (bnfi_marker[k] == NULL)
1334 {
1335 bdr_attr_marker = 1;
1336 break;
1337 }
1338 Array<int> &bdr_marker = *bnfi_marker[k];
1339 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1340 "invalid boundary marker for boundary integrator #"
1341 << k << ", counting from zero");
1342 for (int i = 0; i < bdr_attr_marker.Size(); i++)
1343 {
1344 bdr_attr_marker[i] |= bdr_marker[i];
1345 }
1346 }
1347
1348 for (int i = 0; i < mesh->GetNBE(); i++)
1349 {
1350 const int bdr_attr = mesh->GetBdrAttribute(i);
1351 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1352
1353 T = fes[0]->GetBdrElementTransformation(i);
1354 for (int s = 0; s < fes.Size(); ++s)
1355 {
1356 fe[s] = fes[s]->GetBE(i);
1357 doftrans[s] = fes[s]->GetBdrElementVDofs(i, *(vdofs[s]));
1358 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
1359 if (doftrans[s]) {doftrans[s]->InvTransformPrimal(*el_x[s]); }
1360 }
1361
1362 for (int k = 0; k < bnfi.Size(); k++)
1363 {
1364 if (bnfi_marker[k] &&
1365 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
1366
1367 bnfi[k]->AssembleElementGrad(fe, *T, el_x_const, elmats);
1368
1369 for (int j=0; j<fes.Size(); ++j)
1370 {
1371 for (int l=0; l<fes.Size(); ++l)
1372 {
1373 if (elmats(j,l)->Height() == 0) { continue; }
1374 if (doftrans[j] || doftrans[l])
1375 {
1376 TransformDual(doftrans[j], doftrans[l], *elmats(j,l));
1377 }
1378 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1379 *elmats(j,l), skip_zeros);
1380 }
1381 }
1382 }
1383 }
1384 }
1385
1386 if (fnfi.Size())
1387 {
1389
1390 for (int i = 0; i < mesh->GetNumFaces(); ++i)
1391 {
1392 tr = mesh->GetInteriorFaceTransformations(i);
1393
1394 for (int s=0; s < fes.Size(); ++s)
1395 {
1396 fe[s] = fes[s]->GetFE(tr->Elem1No);
1397 fe2[s] = fes[s]->GetFE(tr->Elem2No);
1398
1399 fes[s]->GetElementVDofs(tr->Elem1No, *vdofs[s]);
1400 fes[s]->GetElementVDofs(tr->Elem2No, *vdofs2[s]);
1401 vdofs[s]->Append(*(vdofs2[s]));
1402
1403 bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
1404 }
1405
1406 for (int k = 0; k < fnfi.Size(); ++k)
1407 {
1408 fnfi[k]->AssembleFaceGrad(fe, fe2, *tr, el_x_const, elmats);
1409 for (int j=0; j<fes.Size(); ++j)
1410 {
1411 for (int l=0; l<fes.Size(); ++l)
1412 {
1413 if (elmats(j,l)->Height() == 0) { continue; }
1414 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1415 *elmats(j,l), skip_zeros);
1416 }
1417 }
1418 }
1419 }
1420 }
1421
1422 if (bfnfi.Size())
1423 {
1425
1426 // Which boundary attributes need to be processed?
1427 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1428 mesh->bdr_attributes.Max() : 0);
1429 bdr_attr_marker = 0;
1430 for (int k = 0; k < bfnfi.Size(); ++k)
1431 {
1432 if (bfnfi_marker[k] == NULL)
1433 {
1434 bdr_attr_marker = 1;
1435 break;
1436 }
1437 Array<int> &bdr_marker = *bfnfi_marker[k];
1438 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1439 "invalid boundary marker for boundary face integrator #"
1440 << k << ", counting from zero");
1441 for (int i = 0; i < bdr_attr_marker.Size(); ++i)
1442 {
1443 bdr_attr_marker[i] |= bdr_marker[i];
1444 }
1445 }
1446
1447 for (int i = 0; i < mesh->GetNBE(); ++i)
1448 {
1449 const int bdr_attr = mesh->GetBdrAttribute(i);
1450 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1451
1452 tr = mesh->GetBdrFaceTransformations(i);
1453 if (tr != NULL)
1454 {
1455 for (int s = 0; s < fes.Size(); ++s)
1456 {
1457 fe[s] = fes[s]->GetFE(tr->Elem1No);
1458 fe2[s] = fe[s];
1459
1460 fes[s]->GetElementVDofs(tr->Elem1No, *vdofs[s]);
1461 bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
1462 }
1463
1464 for (int k = 0; k < bfnfi.Size(); ++k)
1465 {
1466 if (bfnfi_marker[k] &&
1467 (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
1468 bfnfi[k]->AssembleFaceGrad(fe, fe2, *tr, el_x_const, elmats);
1469 for (int l=0; l<fes.Size(); ++l)
1470 {
1471 for (int j=0; j<fes.Size(); ++j)
1472 {
1473 if (elmats(j,l)->Height() == 0) { continue; }
1474 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1475 *elmats(j,l), skip_zeros);
1476 }
1477 }
1478 }
1479 }
1480 }
1481 }
1482
1483 if (!Grads(0,0)->Finalized())
1484 {
1485 for (int i=0; i<fes.Size(); ++i)
1486 {
1487 for (int j=0; j<fes.Size(); ++j)
1488 {
1489 Grads(i,j)->Finalize(skip_zeros);
1490 }
1491 }
1492 }
1493
1494 for (int i=0; i<fes.Size(); ++i)
1495 {
1496 for (int j=0; j<fes.Size(); ++j)
1497 {
1498 delete elmats(i,j);
1499 }
1500 delete vdofs2[i];
1501 delete vdofs[i];
1502 delete el_x[i];
1503 }
1504}
1505
1507{
1508 BlockVector bx(const_cast<Vector&>(x), block_trueOffsets);
1509 const BlockVector &pbx = Prolongate(bx);
1510
1512
1513 Array2D<SparseMatrix *> mGrads(fes.Size(), fes.Size());
1514 mGrads = Grads;
1516 {
1517 for (int s1 = 0; s1 < fes.Size(); ++s1)
1518 {
1519 for (int s2 = 0; s2 < fes.Size(); ++s2)
1520 {
1521 delete cGrads(s1, s2);
1522 cGrads(s1, s2) = RAP(*cP[s1], *Grads(s1, s2), *cP[s2]);
1523 mGrads(s1, s2) = cGrads(s1, s2);
1524 }
1525 }
1526 }
1527
1528 for (int s = 0; s < fes.Size(); ++s)
1529 {
1530 for (int i = 0; i < ess_tdofs[s]->Size(); ++i)
1531 {
1532 for (int j = 0; j < fes.Size(); ++j)
1533 {
1534 if (s == j)
1535 {
1536 mGrads(s, s)->EliminateRowCol((*ess_tdofs[s])[i],
1538 }
1539 else
1540 {
1541 mGrads(s, j)->EliminateRow((*ess_tdofs[s])[i]);
1542 mGrads(j, s)->EliminateCol((*ess_tdofs[s])[i]);
1543 }
1544 }
1545 }
1546 }
1547
1548 delete BlockGrad;
1550 for (int i = 0; i < fes.Size(); ++i)
1551 {
1552 for (int j = 0; j < fes.Size(); ++j)
1553 {
1554 BlockGrad->SetBlock(i, j, mGrads(i, j));
1555 }
1556 }
1557 return *BlockGrad;
1558}
1559
1561{
1562 delete BlockGrad;
1563 for (int i=0; i<fes.Size(); ++i)
1564 {
1565 for (int j=0; j<fes.Size(); ++j)
1566 {
1567 delete Grads(i,j);
1568 delete cGrads(i,j);
1569 }
1570 delete ess_tdofs[i];
1571 }
1572
1573 for (int i = 0; i < dnfi.Size(); ++i)
1574 {
1575 delete dnfi[i];
1576 }
1577
1578 for (int i = 0; i < bnfi.Size(); ++i)
1579 {
1580 delete bnfi[i];
1581 }
1582
1583 for (int i = 0; i < fnfi.Size(); ++i)
1584 {
1585 delete fnfi[i];
1586 }
1587
1588 for (int i = 0; i < bfnfi.Size(); ++i)
1589 {
1590 delete bfnfi[i];
1591 }
1592
1593}
1594
1595}
Dynamic 2D array using row-major layout.
Definition array.hpp:392
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
const T * HostRead() const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), false).
Definition array.hpp:341
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:758
int Size() const
Return the logical size of the array.
Definition array.hpp:147
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
Definition array.cpp:103
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:830
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition array.hpp:935
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:337
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.
void ComputeGradientBlocked(const BlockVector &bx) const
Specialized version of GetGradient() for BlockVector.
real_t GetEnergyBlocked(const BlockVector &bx) const
Specialized version of GetEnergy() for BlockVectors.
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)
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 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:81
void InvTransformPrimal(real_t *v) const
Definition doftrans.cpp:49
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:3881
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:851
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition fespace.hpp:924
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:809
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:658
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition fespace.hpp:900
virtual const Operator * GetProlongationMatrix() const
Definition fespace.hpp:727
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:332
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:3835
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i-th boundary element.
Definition fespace.hpp:933
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:891
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:790
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:822
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:679
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:848
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:348
Abstract class for all finite elements.
Definition fe_base.hpp:244
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:797
Data and methods for unassembled nonlinear forms.
Mesh data type.
Definition mesh.hpp:64
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:290
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition mesh.cpp:6531
int GetAttribute(int i) const
Return the attribute of element i.
Definition mesh.hpp:1389
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition mesh.hpp:1395
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Builds the transformation defining the given boundary face.
Definition mesh.cpp:1123
int GetNBE() const
Returns number of boundary elements.
Definition mesh.hpp:1285
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
See GetFaceElementTransformations().
Definition mesh.cpp:1103
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition mesh.hpp:288
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:29
int GetTrueVSize() const override
Return the number of local vector true dofs.
Definition pfespace.hpp:350
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition pfespace.hpp:388
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:510
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
Definition vector.cpp:679
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:745
int Size() const
Returns the size of the vector.
Definition vector.hpp:226
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition vector.hpp:144
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:558
virtual real_t * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:514
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition vector.cpp:653
void mfem_error(const char *msg)
Definition error.cpp:154
void TransformDual(const DofTransformation *ran_dof_trans, const DofTransformation *dom_dof_trans, DenseMatrix &elmat)
Definition doftrans.cpp:160
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:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
void forall(int N, lambda &&body)
Definition forall.hpp:753