MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
nonlinearform.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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 for (int i = 0; i < ess_tdof_list.Size(); i++)
426 {
427 y(ess_tdof_list[i]) = 0.0;
428 }
429 // y(ess_tdof_list[i]) = x(ess_tdof_list[i]);
430 }
431 // In parallel, the result is in 'py' which is an alias for 'aux2'.
432}
433
435{
436 if (ext)
437 {
438 hGrad.Clear();
439 Operator &grad = ext->GetGradient(Prolongate(x));
440 Operator *Gop;
442 hGrad.Reset(Gop);
443 // In both serial and parallel, when using extension, we return the final
444 // global true-dof gradient with imposed b.c.
445 return *hGrad;
446 }
447
448 const int skip_zeros = 0;
449 Array<int> vdofs;
450 Vector el_x;
451 DenseMatrix elmat;
452 const FiniteElement *fe;
454 DofTransformation *doftrans;
455 Mesh *mesh = fes->GetMesh();
456 const Vector &px = Prolongate(x);
457
458 if (Grad == NULL)
459 {
460 Grad = new SparseMatrix(fes->GetVSize());
461 }
462 else
463 {
464 *Grad = 0.0;
465 }
466
467 if (dnfi.Size())
468 {
469 // Which attributes need to be processed?
470 Array<int> attr_marker(mesh->attributes.Size() ?
471 mesh->attributes.Max() : 0);
472 attr_marker = 0;
473 for (int k = 0; k < dnfi.Size(); k++)
474 {
475 if (dnfi_marker[k] == NULL)
476 {
477 attr_marker = 1;
478 break;
479 }
480 Array<int> &marker = *dnfi_marker[k];
481 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
482 "invalid marker for domain integrator #"
483 << k << ", counting from zero");
484 for (int i = 0; i < attr_marker.Size(); i++)
485 {
486 attr_marker[i] |= marker[i];
487 }
488 }
489
490 for (int i = 0; i < fes->GetNE(); i++)
491 {
492 const int attr = mesh->GetAttribute(i);
493 if (attr_marker[attr-1] == 0) { continue; }
494
495 fe = fes->GetFE(i);
496 doftrans = fes->GetElementVDofs(i, vdofs);
498 px.GetSubVector(vdofs, el_x);
499 if (doftrans) {doftrans->InvTransformPrimal(el_x); }
500 for (int k = 0; k < dnfi.Size(); k++)
501 {
502 if (dnfi_marker[k] &&
503 (*dnfi_marker[k])[attr-1] == 0) { continue; }
504
505 dnfi[k]->AssembleElementGrad(*fe, *T, el_x, elmat);
506 if (doftrans) { doftrans->TransformDual(elmat); }
507 Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
508 // Grad->AddSubMatrix(vdofs, vdofs, elmat, 1);
509 }
510 }
511 }
512
513 if (bnfi.Size())
514 {
515 // Which boundary attributes need to be processed?
516 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
517 mesh->bdr_attributes.Max() : 0);
518 bdr_attr_marker = 0;
519 for (int k = 0; k < bnfi.Size(); k++)
520 {
521 if (bnfi_marker[k] == NULL)
522 {
523 bdr_attr_marker = 1;
524 break;
525 }
526 Array<int> &bdr_marker = *bnfi_marker[k];
527 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
528 "invalid boundary marker for boundary integrator #"
529 << k << ", counting from zero");
530 for (int i = 0; i < bdr_attr_marker.Size(); i++)
531 {
532 bdr_attr_marker[i] |= bdr_marker[i];
533 }
534 }
535
536 for (int i = 0; i < fes->GetNBE(); i++)
537 {
538 const int bdr_attr = mesh->GetBdrAttribute(i);
539 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
540
541 fe = fes->GetBE(i);
542 doftrans = fes->GetBdrElementVDofs(i, vdofs);
544 px.GetSubVector(vdofs, el_x);
545 if (doftrans) {doftrans->InvTransformPrimal(el_x); }
546 for (int k = 0; k < bnfi.Size(); k++)
547 {
548 if (bnfi_marker[k] &&
549 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
550
551 bnfi[k]->AssembleElementGrad(*fe, *T, el_x, elmat);
552 if (doftrans) { doftrans->TransformDual(elmat); }
553 Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
554 }
555 }
556 }
557
558 if (fnfi.Size())
559 {
561 const FiniteElement *fe1, *fe2;
562 Array<int> vdofs2;
563
564 for (int i = 0; i < mesh->GetNumFaces(); i++)
565 {
566 tr = mesh->GetInteriorFaceTransformations(i);
567 if (tr != NULL)
568 {
569 fes->GetElementVDofs(tr->Elem1No, vdofs);
570 fes->GetElementVDofs(tr->Elem2No, vdofs2);
571 vdofs.Append (vdofs2);
572
573 px.GetSubVector(vdofs, el_x);
574
575 fe1 = fes->GetFE(tr->Elem1No);
576 fe2 = fes->GetFE(tr->Elem2No);
577
578 for (int k = 0; k < fnfi.Size(); k++)
579 {
580 fnfi[k]->AssembleFaceGrad(*fe1, *fe2, *tr, el_x, elmat);
581 Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
582 }
583 }
584 }
585 }
586
587 if (bfnfi.Size())
588 {
590 const FiniteElement *fe1, *fe2;
591
592 // Which boundary attributes need to be processed?
593 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
594 mesh->bdr_attributes.Max() : 0);
595 bdr_attr_marker = 0;
596 for (int k = 0; k < bfnfi.Size(); k++)
597 {
598 if (bfnfi_marker[k] == NULL)
599 {
600 bdr_attr_marker = 1;
601 break;
602 }
603 Array<int> &bdr_marker = *bfnfi_marker[k];
604 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
605 "invalid boundary marker for boundary face integrator #"
606 << k << ", counting from zero");
607 for (int i = 0; i < bdr_attr_marker.Size(); i++)
608 {
609 bdr_attr_marker[i] |= bdr_marker[i];
610 }
611 }
612
613 for (int i = 0; i < fes -> GetNBE(); i++)
614 {
615 const int bdr_attr = mesh->GetBdrAttribute(i);
616 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
617
618 tr = mesh->GetBdrFaceTransformations (i);
619 if (tr != NULL)
620 {
621 fes->GetElementVDofs(tr->Elem1No, vdofs);
622 px.GetSubVector(vdofs, el_x);
623
624 fe1 = fes->GetFE(tr->Elem1No);
625 // The fe2 object is really a dummy and not used on the boundaries,
626 // but we can't dereference a NULL pointer, and we don't want to
627 // actually make a fake element.
628 fe2 = fe1;
629 for (int k = 0; k < bfnfi.Size(); k++)
630 {
631 if (bfnfi_marker[k] &&
632 (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
633
634 bfnfi[k]->AssembleFaceGrad(*fe1, *fe2, *tr, el_x, elmat);
635 Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
636 }
637 }
638 }
639 }
640
641 if (!Grad->Finalized())
642 {
643 Grad->Finalize(skip_zeros);
644 }
645
646 SparseMatrix *mGrad = Grad;
647 if (Serial())
648 {
649 if (cP)
650 {
651 delete cGrad;
652 cGrad = RAP(*cP, *Grad, *cP);
653 mGrad = cGrad;
654 }
655 for (int i = 0; i < ess_tdof_list.Size(); i++)
656 {
658 }
659 }
660
661 return *mGrad;
662}
663
665{
666 if (sequence == fes->GetSequence()) { return; }
667
669 delete cGrad; cGrad = NULL;
670 delete Grad; Grad = NULL;
671 hGrad.Clear();
672 ess_tdof_list.SetSize(0); // essential b.c. will need to be set again
674 // Do not modify aux1 and aux2, their size will be set before use.
676 cP = dynamic_cast<const SparseMatrix*>(P);
677
678 if (ext) { ext->Update(); }
679}
680
682{
683 if (ext) { ext->Assemble(); }
684}
685
687{
688 delete cGrad;
689 delete Grad;
690 if (!extern_bfs)
691 {
692 for (int i = 0; i < dnfi.Size(); i++) { delete dnfi[i]; }
693 for (int i = 0; i < bnfi.Size(); i++) { delete bnfi[i]; }
694 for (int i = 0; i < fnfi.Size(); i++) { delete fnfi[i]; }
695 for (int i = 0; i < bfnfi.Size(); i++) { delete bfnfi[i]; }
696 }
697 delete ext;
698}
699
700
702 fes(0), BlockGrad(NULL)
703{
704 height = 0;
705 width = 0;
706}
707
709{
710 delete BlockGrad;
711 BlockGrad = NULL;
712 for (int i=0; i<Grads.NumRows(); ++i)
713 {
714 for (int j=0; j<Grads.NumCols(); ++j)
715 {
716 delete Grads(i,j);
717 delete cGrads(i,j);
718 }
719 }
720 for (int i = 0; i < ess_tdofs.Size(); ++i)
721 {
722 delete ess_tdofs[i];
723 }
724
725 height = 0;
726 width = 0;
727 f.Copy(fes);
728 block_offsets.SetSize(f.Size() + 1);
729 block_trueOffsets.SetSize(f.Size() + 1);
730 block_offsets[0] = 0;
731 block_trueOffsets[0] = 0;
732
733 for (int i=0; i<fes.Size(); ++i)
734 {
735 block_offsets[i+1] = fes[i]->GetVSize();
736 block_trueOffsets[i+1] = fes[i]->GetTrueVSize();
737 }
738
741
742 height = block_trueOffsets[fes.Size()];
743 width = block_trueOffsets[fes.Size()];
744
745 Grads.SetSize(fes.Size(), fes.Size());
746 Grads = NULL;
747
748 cGrads.SetSize(fes.Size(), fes.Size());
749 cGrads = NULL;
750
751 P.SetSize(fes.Size());
752 cP.SetSize(fes.Size());
753 ess_tdofs.SetSize(fes.Size());
754 for (int s = 0; s < fes.Size(); ++s)
755 {
756 // Retrieve prolongation matrix for each FE space
757 P[s] = fes[s]->GetProlongationMatrix();
758 cP[s] = dynamic_cast<const SparseMatrix *>(P[s]);
759
760 // If the P Operator exists and its type is not SparseMatrix, this
761 // indicates the Operator is part of parallel run.
762 if (P[s] && !cP[s])
763 {
764 is_serial = false;
765 }
766
767 // If the P Operator exists and its type is SparseMatrix, this indicates
768 // the Operator is serial but needs prolongation on assembly.
769 if (cP[s])
770 {
771 needs_prolongation = true;
772 }
773
774 ess_tdofs[s] = new Array<int>;
775 }
776}
777
783
785 const Array<Array<int> *> &bdr_attr_is_ess, Array<Vector *> &rhs)
786{
787 for (int s = 0; s < fes.Size(); ++s)
788 {
789 ess_tdofs[s]->SetSize(ess_tdofs.Size());
790
791 fes[s]->GetEssentialTrueDofs(*bdr_attr_is_ess[s], *ess_tdofs[s]);
792
793 if (rhs[s])
794 {
795 rhs[s]->SetSubVector(*ess_tdofs[s], 0.0);
796 }
797 }
798}
799
801{
802 Array<Array<int> *> vdofs(fes.Size());
803 Array<Vector *> el_x(fes.Size());
804 Array<const Vector *> el_x_const(fes.Size());
807 DofTransformation *doftrans;
808 Mesh *mesh = fes[0]->GetMesh();
809 real_t energy = 0.0;
810
811 for (int i=0; i<fes.Size(); ++i)
812 {
813 el_x_const[i] = el_x[i] = new Vector();
814 vdofs[i] = new Array<int>;
815 }
816
817 if (dnfi.Size())
818 {
819 // Which attributes need to be processed?
820 Array<int> attr_marker(mesh->attributes.Size() ?
821 mesh->attributes.Max() : 0);
822 attr_marker = 0;
823 for (int k = 0; k < dnfi.Size(); k++)
824 {
825 if (dnfi_marker[k] == NULL)
826 {
827 attr_marker = 1;
828 break;
829 }
830 Array<int> &marker = *dnfi_marker[k];
831 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
832 "invalid marker for domain integrator #"
833 << k << ", counting from zero");
834 for (int i = 0; i < attr_marker.Size(); i++)
835 {
836 attr_marker[i] |= marker[i];
837 }
838 }
839
840 for (int i = 0; i < fes[0]->GetNE(); ++i)
841 {
842 const int attr = mesh->GetAttribute(i);
843 if (attr_marker[attr-1] == 0) { continue; }
844
845 T = fes[0]->GetElementTransformation(i);
846 for (int s=0; s<fes.Size(); ++s)
847 {
848 fe[s] = fes[s]->GetFE(i);
849 doftrans = fes[s]->GetElementVDofs(i, *vdofs[s]);
850 bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
851 if (doftrans) {doftrans->InvTransformPrimal(*el_x[s]); }
852 }
853
854 for (int k = 0; k < dnfi.Size(); ++k)
855 {
856 if (dnfi_marker[k] &&
857 (*dnfi_marker[k])[attr-1] == 0) { continue; }
858
859 energy += dnfi[k]->GetElementEnergy(fe, *T, el_x_const);
860 }
861 }
862 }
863
864 if (bnfi.Size())
865 {
866 // Which boundary attributes need to be processed?
867 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
868 mesh->bdr_attributes.Max() : 0);
869 bdr_attr_marker = 0;
870 for (int k = 0; k < bnfi.Size(); k++)
871 {
872 if (bnfi_marker[k] == NULL)
873 {
874 bdr_attr_marker = 1;
875 break;
876 }
877 Array<int> &bdr_marker = *bnfi_marker[k];
878 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
879 "invalid boundary marker for boundary integrator #"
880 << k << ", counting from zero");
881 for (int i = 0; i < bdr_attr_marker.Size(); i++)
882 {
883 bdr_attr_marker[i] |= bdr_marker[i];
884 }
885 }
886
887 for (int i = 0; i < mesh->GetNBE(); i++)
888 {
889 const int bdr_attr = mesh->GetBdrAttribute(i);
890 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
891
892 T = fes[0]->GetBdrElementTransformation(i);
893 for (int s = 0; s < fes.Size(); ++s)
894 {
895 fe[s] = fes[s]->GetBE(i);
896 doftrans = fes[s]->GetBdrElementVDofs(i, *(vdofs[s]));
897 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
898 if (doftrans) {doftrans->InvTransformPrimal(*el_x[s]); }
899 }
900
901 for (int k = 0; k < bnfi.Size(); k++)
902 {
903 if (bnfi_marker[k] &&
904 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
905
906 energy += bnfi[k]->GetElementEnergy(fe, *T, el_x_const);
907 }
908 }
909 }
910
911 // free the allocated memory
912 for (int i = 0; i < fes.Size(); ++i)
913 {
914 delete el_x[i];
915 delete vdofs[i];
916 }
917
918 if (fnfi.Size())
919 {
920 MFEM_ABORT("TODO: add energy contribution from interior face terms");
921 }
922
923 if (bfnfi.Size())
924 {
925 MFEM_ABORT("TODO: add energy contribution from boundary face terms");
926 }
927
928 return energy;
929}
930
932{
933 xs.Update(const_cast<Vector&>(x), block_offsets);
934 return GetEnergyBlocked(xs);
935}
936
938 BlockVector &by) const
939{
940 Array<Array<int> *>vdofs(fes.Size());
941 Array<Array<int> *>vdofs2(fes.Size());
942 Array<Vector *> el_x(fes.Size());
943 Array<const Vector *> el_x_const(fes.Size());
944 Array<Vector *> el_y(fes.Size());
948 Array<DofTransformation *> doftrans(fes.Size()); doftrans = nullptr;
949 Mesh *mesh = fes[0]->GetMesh();
950
951 by.UseDevice(true);
952 by = 0.0;
953 by.SyncToBlocks();
954 for (int s=0; s<fes.Size(); ++s)
955 {
956 el_x_const[s] = el_x[s] = new Vector();
957 el_y[s] = new Vector();
958 vdofs[s] = new Array<int>;
959 vdofs2[s] = new Array<int>;
960 }
961
962 if (dnfi.Size())
963 {
964 // Which attributes need to be processed?
965 Array<int> attr_marker(mesh->attributes.Size() ?
966 mesh->attributes.Max() : 0);
967 attr_marker = 0;
968 for (int k = 0; k < dnfi.Size(); k++)
969 {
970 if (dnfi_marker[k] == NULL)
971 {
972 attr_marker = 1;
973 break;
974 }
975 Array<int> &marker = *dnfi_marker[k];
976 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
977 "invalid marker for domain integrator #"
978 << k << ", counting from zero");
979 for (int i = 0; i < attr_marker.Size(); i++)
980 {
981 attr_marker[i] |= marker[i];
982 }
983 }
984
985 for (int i = 0; i < fes[0]->GetNE(); ++i)
986 {
987 const int attr = mesh->GetAttribute(i);
988 if (attr_marker[attr-1] == 0) { continue; }
989
990 T = fes[0]->GetElementTransformation(i);
991 for (int s = 0; s < fes.Size(); ++s)
992 {
993 doftrans[s] = fes[s]->GetElementVDofs(i, *(vdofs[s]));
994 fe[s] = fes[s]->GetFE(i);
995 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
996 if (doftrans[s]) {doftrans[s]->InvTransformPrimal(*el_x[s]); }
997 }
998
999 for (int k = 0; k < dnfi.Size(); ++k)
1000 {
1001 if (dnfi_marker[k] &&
1002 (*dnfi_marker[k])[attr-1] == 0) { continue; }
1003
1004 dnfi[k]->AssembleElementVector(fe, *T,
1005 el_x_const, el_y);
1006
1007 for (int s=0; s<fes.Size(); ++s)
1008 {
1009 if (el_y[s]->Size() == 0) { continue; }
1010 if (doftrans[s]) {doftrans[s]->TransformDual(*el_y[s]); }
1011 by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
1012 }
1013 }
1014 }
1015 }
1016
1017 if (bnfi.Size())
1018 {
1019 // Which boundary attributes need to be processed?
1020 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1021 mesh->bdr_attributes.Max() : 0);
1022 bdr_attr_marker = 0;
1023 for (int k = 0; k < bnfi.Size(); k++)
1024 {
1025 if (bnfi_marker[k] == NULL)
1026 {
1027 bdr_attr_marker = 1;
1028 break;
1029 }
1030 Array<int> &bdr_marker = *bnfi_marker[k];
1031 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1032 "invalid boundary marker for boundary integrator #"
1033 << k << ", counting from zero");
1034 for (int i = 0; i < bdr_attr_marker.Size(); i++)
1035 {
1036 bdr_attr_marker[i] |= bdr_marker[i];
1037 }
1038 }
1039
1040 for (int i = 0; i < mesh->GetNBE(); i++)
1041 {
1042 const int bdr_attr = mesh->GetBdrAttribute(i);
1043 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1044
1045 T = fes[0]->GetBdrElementTransformation(i);
1046 for (int s = 0; s < fes.Size(); ++s)
1047 {
1048 doftrans[s] = fes[s]->GetBdrElementVDofs(i, *(vdofs[s]));
1049 fe[s] = fes[s]->GetBE(i);
1050 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
1051 if (doftrans[s]) {doftrans[s]->InvTransformPrimal(*el_x[s]); }
1052 }
1053
1054 for (int k = 0; k < bnfi.Size(); k++)
1055 {
1056 if (bnfi_marker[k] &&
1057 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
1058
1059 bnfi[k]->AssembleElementVector(fe, *T, el_x_const, el_y);
1060
1061 for (int s=0; s<fes.Size(); ++s)
1062 {
1063 if (el_y[s]->Size() == 0) { continue; }
1064 if (doftrans[s]) {doftrans[s]->TransformDual(*el_y[s]); }
1065 by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
1066 }
1067 }
1068 }
1069 }
1070
1071
1072 if (fnfi.Size())
1073 {
1075
1076 for (int i = 0; i < mesh->GetNumFaces(); ++i)
1077 {
1078 tr = mesh->GetInteriorFaceTransformations(i);
1079 if (tr != NULL)
1080 {
1081 for (int s=0; s<fes.Size(); ++s)
1082 {
1083 fe[s] = fes[s]->GetFE(tr->Elem1No);
1084 fe2[s] = fes[s]->GetFE(tr->Elem2No);
1085
1086 fes[s]->GetElementVDofs(tr->Elem1No, *(vdofs[s]));
1087 fes[s]->GetElementVDofs(tr->Elem2No, *(vdofs2[s]));
1088
1089 vdofs[s]->Append(*(vdofs2[s]));
1090
1091 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
1092 }
1093
1094 for (int k = 0; k < fnfi.Size(); ++k)
1095 {
1096
1097 fnfi[k]->AssembleFaceVector(fe, fe2, *tr, el_x_const, el_y);
1098
1099 for (int s=0; s<fes.Size(); ++s)
1100 {
1101 if (el_y[s]->Size() == 0) { continue; }
1102 by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
1103 }
1104 }
1105 }
1106 }
1107 }
1108
1109 if (bfnfi.Size())
1110 {
1112
1113 // Which boundary attributes need to be processed?
1114 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1115 mesh->bdr_attributes.Max() : 0);
1116 bdr_attr_marker = 0;
1117 for (int k = 0; k < bfnfi.Size(); ++k)
1118 {
1119 if (bfnfi_marker[k] == NULL)
1120 {
1121 bdr_attr_marker = 1;
1122 break;
1123 }
1124 Array<int> &bdr_marker = *bfnfi_marker[k];
1125 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1126 "invalid boundary marker for boundary face integrator #"
1127 << k << ", counting from zero");
1128 for (int i = 0; i < bdr_attr_marker.Size(); ++i)
1129 {
1130 bdr_attr_marker[i] |= bdr_marker[i];
1131 }
1132 }
1133
1134 for (int i = 0; i < mesh->GetNBE(); ++i)
1135 {
1136 const int bdr_attr = mesh->GetBdrAttribute(i);
1137 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1138
1139 tr = mesh->GetBdrFaceTransformations(i);
1140 if (tr != NULL)
1141 {
1142 for (int s=0; s<fes.Size(); ++s)
1143 {
1144 fe[s] = fes[s]->GetFE(tr->Elem1No);
1145 fe2[s] = fes[s]->GetFE(tr->Elem1No);
1146
1147 fes[s]->GetElementVDofs(tr->Elem1No, *(vdofs[s]));
1148 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
1149 }
1150
1151 for (int k = 0; k < bfnfi.Size(); ++k)
1152 {
1153 if (bfnfi_marker[k] &&
1154 (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
1155
1156 bfnfi[k]->AssembleFaceVector(fe, fe2, *tr, el_x_const, el_y);
1157
1158 for (int s=0; s<fes.Size(); ++s)
1159 {
1160 if (el_y[s]->Size() == 0) { continue; }
1161 by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
1162 }
1163 }
1164 }
1165 }
1166 }
1167
1168 for (int s=0; s<fes.Size(); ++s)
1169 {
1170 delete vdofs2[s];
1171 delete vdofs[s];
1172 delete el_y[s];
1173 delete el_x[s];
1174 }
1175
1176 by.SyncFromBlocks();
1177}
1178
1180{
1181 MFEM_VERIFY(bx.Size() == Width(), "invalid input BlockVector size");
1182
1184 {
1186 for (int s = 0; s < fes.Size(); s++)
1187 {
1188 P[s]->Mult(bx.GetBlock(s), aux1.GetBlock(s));
1189 }
1190 return aux1;
1191 }
1192 return bx;
1193}
1194
1196{
1197 BlockVector bx(const_cast<Vector&>(x), block_trueOffsets);
1199
1200 const BlockVector &pbx = Prolongate(bx);
1202 {
1204 }
1205 BlockVector &pby = needs_prolongation ? aux2 : by;
1206
1207 xs.Update(const_cast<BlockVector&>(pbx), block_offsets);
1208 ys.Update(pby, block_offsets);
1209 MultBlocked(xs, ys);
1210
1211 for (int s = 0; s < fes.Size(); s++)
1212 {
1213 if (cP[s])
1214 {
1215 cP[s]->MultTranspose(pby.GetBlock(s), by.GetBlock(s));
1216 }
1217 by.GetBlock(s).SetSubVector(*ess_tdofs[s], 0.0);
1218 }
1219}
1220
1222{
1223 const int skip_zeros = 0;
1224 Array<Array<int> *> vdofs(fes.Size());
1225 Array<Array<int> *> vdofs2(fes.Size());
1226 Array<Vector *> el_x(fes.Size());
1227 Array<const Vector *> el_x_const(fes.Size());
1228 Array2D<DenseMatrix *> elmats(fes.Size(), fes.Size());
1232 Array<DofTransformation *> doftrans(fes.Size()); doftrans = nullptr;
1233 Mesh *mesh = fes[0]->GetMesh();
1234
1235 for (int i=0; i<fes.Size(); ++i)
1236 {
1237 el_x_const[i] = el_x[i] = new Vector();
1238 vdofs[i] = new Array<int>;
1239 vdofs2[i] = new Array<int>;
1240 for (int j=0; j<fes.Size(); ++j)
1241 {
1242 elmats(i,j) = new DenseMatrix();
1243 }
1244 }
1245
1246 for (int i=0; i<fes.Size(); ++i)
1247 {
1248 for (int j=0; j<fes.Size(); ++j)
1249 {
1250 if (Grads(i,j) != NULL)
1251 {
1252 *Grads(i,j) = 0.0;
1253 }
1254 else
1255 {
1256 Grads(i,j) = new SparseMatrix(fes[i]->GetVSize(),
1257 fes[j]->GetVSize());
1258 }
1259 }
1260 }
1261
1262 if (dnfi.Size())
1263 {
1264 // Which attributes need to be processed?
1265 Array<int> attr_marker(mesh->attributes.Size() ?
1266 mesh->attributes.Max() : 0);
1267 attr_marker = 0;
1268 for (int k = 0; k < dnfi.Size(); k++)
1269 {
1270 if (dnfi_marker[k] == NULL)
1271 {
1272 attr_marker = 1;
1273 break;
1274 }
1275 Array<int> &marker = *dnfi_marker[k];
1276 MFEM_ASSERT(marker.Size() == attr_marker.Size(),
1277 "invalid marker for domain integrator #"
1278 << k << ", counting from zero");
1279 for (int i = 0; i < attr_marker.Size(); i++)
1280 {
1281 attr_marker[i] |= marker[i];
1282 }
1283 }
1284
1285 for (int i = 0; i < fes[0]->GetNE(); ++i)
1286 {
1287 const int attr = mesh->GetAttribute(i);
1288 if (attr_marker[attr-1] == 0) { continue; }
1289
1290 T = fes[0]->GetElementTransformation(i);
1291 for (int s = 0; s < fes.Size(); ++s)
1292 {
1293 fe[s] = fes[s]->GetFE(i);
1294 doftrans[s] = fes[s]->GetElementVDofs(i, *vdofs[s]);
1295 bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
1296 if (doftrans[s]) {doftrans[s]->InvTransformPrimal(*el_x[s]); }
1297 }
1298
1299 for (int k = 0; k < dnfi.Size(); ++k)
1300 {
1301 if (dnfi_marker[k] &&
1302 (*dnfi_marker[k])[attr-1] == 0) { continue; }
1303
1304 dnfi[k]->AssembleElementGrad(fe, *T, el_x_const, elmats);
1305
1306 for (int j=0; j<fes.Size(); ++j)
1307 {
1308 for (int l=0; l<fes.Size(); ++l)
1309 {
1310 if (elmats(j,l)->Height() == 0) { continue; }
1311 if (doftrans[j] || doftrans[l])
1312 {
1313 TransformDual(doftrans[j], doftrans[l], *elmats(j,l));
1314 }
1315 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1316 *elmats(j,l), skip_zeros);
1317 }
1318 }
1319 }
1320 }
1321 }
1322
1323 if (bnfi.Size())
1324 {
1325 // Which boundary attributes need to be processed?
1326 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1327 mesh->bdr_attributes.Max() : 0);
1328 bdr_attr_marker = 0;
1329 for (int k = 0; k < bnfi.Size(); k++)
1330 {
1331 if (bnfi_marker[k] == NULL)
1332 {
1333 bdr_attr_marker = 1;
1334 break;
1335 }
1336 Array<int> &bdr_marker = *bnfi_marker[k];
1337 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1338 "invalid boundary marker for boundary integrator #"
1339 << k << ", counting from zero");
1340 for (int i = 0; i < bdr_attr_marker.Size(); i++)
1341 {
1342 bdr_attr_marker[i] |= bdr_marker[i];
1343 }
1344 }
1345
1346 for (int i = 0; i < mesh->GetNBE(); i++)
1347 {
1348 const int bdr_attr = mesh->GetBdrAttribute(i);
1349 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1350
1351 T = fes[0]->GetBdrElementTransformation(i);
1352 for (int s = 0; s < fes.Size(); ++s)
1353 {
1354 fe[s] = fes[s]->GetBE(i);
1355 doftrans[s] = fes[s]->GetBdrElementVDofs(i, *(vdofs[s]));
1356 bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
1357 if (doftrans[s]) {doftrans[s]->InvTransformPrimal(*el_x[s]); }
1358 }
1359
1360 for (int k = 0; k < bnfi.Size(); k++)
1361 {
1362 if (bnfi_marker[k] &&
1363 (*bnfi_marker[k])[bdr_attr-1] == 0) { continue; }
1364
1365 bnfi[k]->AssembleElementGrad(fe, *T, el_x_const, elmats);
1366
1367 for (int j=0; j<fes.Size(); ++j)
1368 {
1369 for (int l=0; l<fes.Size(); ++l)
1370 {
1371 if (elmats(j,l)->Height() == 0) { continue; }
1372 if (doftrans[j] || doftrans[l])
1373 {
1374 TransformDual(doftrans[j], doftrans[l], *elmats(j,l));
1375 }
1376 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1377 *elmats(j,l), skip_zeros);
1378 }
1379 }
1380 }
1381 }
1382 }
1383
1384 if (fnfi.Size())
1385 {
1387
1388 for (int i = 0; i < mesh->GetNumFaces(); ++i)
1389 {
1390 tr = mesh->GetInteriorFaceTransformations(i);
1391
1392 for (int s=0; s < fes.Size(); ++s)
1393 {
1394 fe[s] = fes[s]->GetFE(tr->Elem1No);
1395 fe2[s] = fes[s]->GetFE(tr->Elem2No);
1396
1397 fes[s]->GetElementVDofs(tr->Elem1No, *vdofs[s]);
1398 fes[s]->GetElementVDofs(tr->Elem2No, *vdofs2[s]);
1399 vdofs[s]->Append(*(vdofs2[s]));
1400
1401 bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
1402 }
1403
1404 for (int k = 0; k < fnfi.Size(); ++k)
1405 {
1406 fnfi[k]->AssembleFaceGrad(fe, fe2, *tr, el_x_const, elmats);
1407 for (int j=0; j<fes.Size(); ++j)
1408 {
1409 for (int l=0; l<fes.Size(); ++l)
1410 {
1411 if (elmats(j,l)->Height() == 0) { continue; }
1412 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1413 *elmats(j,l), skip_zeros);
1414 }
1415 }
1416 }
1417 }
1418 }
1419
1420 if (bfnfi.Size())
1421 {
1423
1424 // Which boundary attributes need to be processed?
1425 Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1426 mesh->bdr_attributes.Max() : 0);
1427 bdr_attr_marker = 0;
1428 for (int k = 0; k < bfnfi.Size(); ++k)
1429 {
1430 if (bfnfi_marker[k] == NULL)
1431 {
1432 bdr_attr_marker = 1;
1433 break;
1434 }
1435 Array<int> &bdr_marker = *bfnfi_marker[k];
1436 MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1437 "invalid boundary marker for boundary face integrator #"
1438 << k << ", counting from zero");
1439 for (int i = 0; i < bdr_attr_marker.Size(); ++i)
1440 {
1441 bdr_attr_marker[i] |= bdr_marker[i];
1442 }
1443 }
1444
1445 for (int i = 0; i < mesh->GetNBE(); ++i)
1446 {
1447 const int bdr_attr = mesh->GetBdrAttribute(i);
1448 if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1449
1450 tr = mesh->GetBdrFaceTransformations(i);
1451 if (tr != NULL)
1452 {
1453 for (int s = 0; s < fes.Size(); ++s)
1454 {
1455 fe[s] = fes[s]->GetFE(tr->Elem1No);
1456 fe2[s] = fe[s];
1457
1458 fes[s]->GetElementVDofs(tr->Elem1No, *vdofs[s]);
1459 bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
1460 }
1461
1462 for (int k = 0; k < bfnfi.Size(); ++k)
1463 {
1464 if (bfnfi_marker[k] &&
1465 (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
1466 bfnfi[k]->AssembleFaceGrad(fe, fe2, *tr, el_x_const, elmats);
1467 for (int l=0; l<fes.Size(); ++l)
1468 {
1469 for (int j=0; j<fes.Size(); ++j)
1470 {
1471 if (elmats(j,l)->Height() == 0) { continue; }
1472 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1473 *elmats(j,l), skip_zeros);
1474 }
1475 }
1476 }
1477 }
1478 }
1479 }
1480
1481 if (!Grads(0,0)->Finalized())
1482 {
1483 for (int i=0; i<fes.Size(); ++i)
1484 {
1485 for (int j=0; j<fes.Size(); ++j)
1486 {
1487 Grads(i,j)->Finalize(skip_zeros);
1488 }
1489 }
1490 }
1491
1492 for (int i=0; i<fes.Size(); ++i)
1493 {
1494 for (int j=0; j<fes.Size(); ++j)
1495 {
1496 delete elmats(i,j);
1497 }
1498 delete vdofs2[i];
1499 delete vdofs[i];
1500 delete el_x[i];
1501 }
1502}
1503
1505{
1506 BlockVector bx(const_cast<Vector&>(x), block_trueOffsets);
1507 const BlockVector &pbx = Prolongate(bx);
1508
1510
1511 Array2D<SparseMatrix *> mGrads(fes.Size(), fes.Size());
1512 mGrads = Grads;
1514 {
1515 for (int s1 = 0; s1 < fes.Size(); ++s1)
1516 {
1517 for (int s2 = 0; s2 < fes.Size(); ++s2)
1518 {
1519 delete cGrads(s1, s2);
1520 cGrads(s1, s2) = RAP(*cP[s1], *Grads(s1, s2), *cP[s2]);
1521 mGrads(s1, s2) = cGrads(s1, s2);
1522 }
1523 }
1524 }
1525
1526 for (int s = 0; s < fes.Size(); ++s)
1527 {
1528 for (int i = 0; i < ess_tdofs[s]->Size(); ++i)
1529 {
1530 for (int j = 0; j < fes.Size(); ++j)
1531 {
1532 if (s == j)
1533 {
1534 mGrads(s, s)->EliminateRowCol((*ess_tdofs[s])[i],
1536 }
1537 else
1538 {
1539 mGrads(s, j)->EliminateRow((*ess_tdofs[s])[i]);
1540 mGrads(j, s)->EliminateCol((*ess_tdofs[s])[i]);
1541 }
1542 }
1543 }
1544 }
1545
1546 delete BlockGrad;
1548 for (int i = 0; i < fes.Size(); ++i)
1549 {
1550 for (int j = 0; j < fes.Size(); ++j)
1551 {
1552 BlockGrad->SetBlock(i, j, mGrads(i, j));
1553 }
1554 }
1555 return *BlockGrad;
1556}
1557
1559{
1560 delete BlockGrad;
1561 for (int i=0; i<fes.Size(); ++i)
1562 {
1563 for (int j=0; j<fes.Size(); ++j)
1564 {
1565 delete Grads(i,j);
1566 delete cGrads(i,j);
1567 }
1568 delete ess_tdofs[i];
1569 }
1570
1571 for (int i = 0; i < dnfi.Size(); ++i)
1572 {
1573 delete dnfi[i];
1574 }
1575
1576 for (int i = 0; i < bnfi.Size(); ++i)
1577 {
1578 delete bnfi[i];
1579 }
1580
1581 for (int i = 0; i < fnfi.Size(); ++i)
1582 {
1583 delete fnfi[i];
1584 }
1585
1586 for (int i = 0; i < bfnfi.Size(); ++i)
1587 {
1588 delete bfnfi[i];
1589 }
1590
1591}
1592
1593}
Dynamic 2D array using row-major layout.
Definition array.hpp:372
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
int Size() const
Return the logical size of the array.
Definition array.hpp:144
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:769
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition array.hpp:874
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:317
Array< Array< int > * > dnfi_marker
Array< Array< int > * > bnfi_marker
Array< BlockNonlinearFormIntegrator * > dnfi
Set of Domain Integrators to be assembled (added).
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
virtual Operator & GetGradient(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().
virtual void Mult(const Vector &x, Vector &y) const
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:484
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th boundary fac...
Definition fespace.cpp:3204
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:716
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition fespace.hpp:773
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:667
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:588
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition fespace.hpp:749
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition fespace.hpp:597
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:287
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:3168
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i-th boundary element.
Definition fespace.hpp:782
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:740
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:648
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:680
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:559
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:713
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:303
Abstract class for all finite elements.
Definition fe_base.hpp:239
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:777
Data and methods for unassembled nonlinear forms.
Mesh data type.
Definition mesh.hpp:56
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:282
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition mesh.cpp:6250
int GetAttribute(int i) const
Return the attribute of element i.
Definition mesh.hpp:1333
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition mesh.hpp:1339
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Builds the transformation defining the given boundary face.
Definition mesh.cpp:1099
int GetNBE() const
Returns number of boundary elements.
Definition mesh.hpp:1229
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
See GetFaceElementTransformations().
Definition mesh.cpp:1079
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition mesh.hpp:280
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.
virtual void Mult(const Vector &x, Vector &y) const
Evaluate the action of the NonlinearForm.
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.
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 Operator & GetGradient(const Vector &x) const
Compute the gradient Operator of the NonlinearForm corresponding to the state x.
virtual ~NonlinearForm()
Destroy the NonlinearForm including the owned NonlinearFormIntegrators and gradient Operator.
SparseMatrix * cGrad
Array< Array< int > * > bnfi_marker
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:289
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition pfespace.hpp:327
Data type sparse matrix.
Definition sparsemat.hpp:51
virtual void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
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.
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)
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR.
Vector data type.
Definition vector.hpp:80
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:490
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
Definition vector.cpp:604
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:670
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition vector.hpp:136
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition vector.cpp:578
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:754
RefCoord s[3]