Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
BCBuilder_def.hpp
1#ifndef BCBuilder_def_hpp
2#define BCBuilder_def_hpp
3
4#include "feddlib/core/Utils/FEDDUtils.hpp"
5#include <Teuchos_RCPDecl.hpp>
6#include <Teuchos_ScalarTraitsDecl.hpp>
7
8
9void dummyFuncBC(double* x, double* res, double t, const double* parameters)
10{
11 return;
12}
13
22
23namespace FEDD {
24template<class SC,class LO,class GO,class NO>
25BCBuilder<SC,LO,GO,NO>::BCBuilder():
26 vecBC_func_(),
27 vecFlag_(),
28 vecBlockID_(),
29 vecDomain_(),
30 vecBCType_(),
31 vecDofs_(),
32 vecBC_Parameters_(),
33 vecExternalSol_(0),
34 resultPtr_(),
35 pointPtr_()
36#ifdef BCBuilder_TIMER
37,SetSystemRowTimer_(Teuchos::TimeMonitor::getNewCounter("BCBuilder: SetSystemRow")),
38BlockRowHasDirichletTimer_(Teuchos::TimeMonitor::getNewCounter("BCBuilder: SetSystemRow: BlockHasDiri")),
39FindFlagTimer_(Teuchos::TimeMonitor::getNewCounter("BCBuilder: SetSystemRow: BlockHasDiri: FindFlag"))
40#endif
41{
42
43}
44
45template<class SC,class LO,class GO,class NO>
46void BCBuilder<SC,LO,GO,NO>::addBC(BC_func_Type funcBC, int flag, int block, const DomainPtr_Type &domain, std::string type, int dofs){
47
48 vecBC_func_.push_back(funcBC);
49 vecFlag_.push_back(flag);
50 vecBlockID_.push_back(block);
51 vecDomain_.push_back(domain);
52 vecBCType_.push_back(type);
53 vecDofs_.push_back(dofs);
54 vec_dbl_Type dummy(1,0.);
55 vecBC_Parameters_.push_back(dummy);
56 MultiVectorConstPtr_Type dummyMV;
57 vecExternalSol_.push_back(dummyMV);
58 vecFlowRateBool_.push_back(false);
59 vecBC_func_flowRate_.push_back(dummyFuncBC);
60
61
62}
63
64template<class SC,class LO,class GO,class NO>
65void BCBuilder<SC,LO,GO,NO>::addBC(BC_func_Type funcBC, int flag, int block, const DomainPtr_Type &domain, std::string type, int dofs, vec_dbl_Type &parameter_vec){
66
67 vecBC_func_.push_back(funcBC);
68 vecFlag_.push_back(flag);
69 vecBlockID_.push_back(block);
70 vecDomain_.push_back(domain);
71 vecBCType_.push_back(type);
72 vecDofs_.push_back(dofs);
73 vecBC_Parameters_.push_back(parameter_vec);
74 MultiVectorConstPtr_Type dummyMV;
75 vecExternalSol_.push_back(dummyMV);
76 vecFlowRateBool_.push_back(false);
77 vecBC_func_flowRate_.push_back(dummyFuncBC);
78
79}
80
81template<class SC,class LO,class GO,class NO>
82void BCBuilder<SC,LO,GO,NO>::addBC(BC_func_Type funcBC, int flag, int block, const DomainPtr_Type &domain, std::string type, int dofs, vec_dbl_Type &parameter_vec, MultiVectorConstPtr_Type& externalSol){
83
84 vecBC_func_.push_back(funcBC);
85 vecFlag_.push_back(flag);
86 vecBlockID_.push_back(block);
87 vecDomain_.push_back(domain);
88 vecBCType_.push_back(type);
89 vecDofs_.push_back(dofs);
90 vecBC_Parameters_.push_back(parameter_vec);
91 vecExternalSol_.push_back(externalSol);
92 vecFlowRateBool_.push_back(false);
93 vecBC_func_flowRate_.push_back(dummyFuncBC);
94}
95
96template<class SC,class LO,class GO,class NO>
97void BCBuilder<SC,LO,GO,NO>::addBC(BC_func_Type funcBC, int flag, int block, const DomainPtr_Type &domain,
98 std::string type, int dofs, vec_dbl_Type &parameter_vec, MultiVectorConstPtr_Type& externalSol,
99 bool determineFlowRate,BC_func_Type funcBC_flowRate){
101 vecBC_func_.push_back(funcBC);
102 vecFlag_.push_back(flag);
103 vecBlockID_.push_back(block);
104 vecDomain_.push_back(domain);
105 vecBCType_.push_back(type);
106 vecDofs_.push_back(dofs);
107 vecBC_Parameters_.push_back(parameter_vec);
108 vecExternalSol_.push_back(externalSol);
109 vecFlowRateBool_.push_back(true);
110 vecBC_func_flowRate_.push_back(funcBC_flowRate);
111}
112
113template<class SC,class LO,class GO,class NO>
114void BCBuilder<SC,LO,GO,NO>::set(const BlockMatrixPtr_Type &blockMatrix, const BlockMultiVectorPtr_Type &blockMV, double t) const{
115
116 setSystem(blockMatrix);
117
118 setRHS(blockMV, t);
119
120}
121
122template<class SC,class LO,class GO,class NO>
123void BCBuilder<SC,LO,GO,NO>::setRHS(const BlockMultiVectorPtr_Type &blockMV, double t) const{
124 vec_dbl_Type result;
125 vec_dbl_Type point;
126
127 TEUCHOS_TEST_FOR_EXCEPTION( blockMV->getNumVectors()>1, std::runtime_error, "BCBuilder::setRHS() only for getNumVectors == 1.");
128
129 int loc = 0;
130 for (int block = 0; block < blockMV->size(); block++) { // blocks of RHS vector
131 if(blockHasDirichletBC(block, loc)){
132 if (vecFlowRateBool_[loc]) { // we use the external vector here with flowrate
133 this->determineVelocityForFlowrate(loc,t);
134 }
135 vec_int_ptr_Type bcFlags = vecDomain_.at(loc)->getBCFlagUnique();
136 vec2D_dbl_ptr_Type vecPoints = vecDomain_.at(loc)->getPointsUnique();
137 int dim = vecDomain_.at(loc)->getDimension();
138 int dofs = vecDofs_.at(loc);
139
140 result.resize(dofs,0.);
141 point.resize(dim,0.);
142
143 resultPtr_ = Teuchos::rcp(new vec_dbl_Type(dofs,0.));
144 pointPtr_ = Teuchos::rcp(new vec_dbl_Type(dim,0.));
145
146 Teuchos::ArrayRCP<SC> valuesRHS = blockMV->getBlock(block)->getDataNonConst(0);
147 for (int i = 0 ; i < bcFlags->size(); i++) {
148 int flag = bcFlags->at(i);
149 if(findFlag(flag, block, loc)){ // loc changed to postion corresponding to bcFlag
150 if (!vecBCType_.at(loc).compare("Dirichlet") ||
151 !vecBCType_.at(loc).compare("Dirichlet_X") ||
152 !vecBCType_.at(loc).compare("Dirichlet_Y") ||
153 !vecBCType_.at(loc).compare("Dirichlet_Z")||
154 !vecBCType_.at(loc).compare("Dirichlet_X_Y") ||
155 !vecBCType_.at(loc).compare("Dirichlet_X_Z") ||
156 !vecBCType_.at(loc).compare("Dirichlet_Y_Z")) {
157 if (!vecExternalSol_[loc].is_null()) { // we use the external vector here
158 std::string type = "standard";
159 this->setDirichletBoundaryFromExternal(valuesRHS, i, loc, t, type);
160 }
161 else {
162 for (int d=0; d < dim; d++) {
163 point.at(d) = vecPoints->at(i).at(d);
164 }
165 for (int d=0; d < dofs; d++) {
166 result.at(d) = vecPoints->at(i).at(d);
167 }
168 vecBC_func_.at(loc)(&point.at(0), &result.at(0), t, &(vecBC_Parameters_.at(loc).at(0)));
170 for (UN j=0; j<blockMV->getNumVectors(); j++) {
171 Teuchos::ArrayRCP<SC> values = blockMV->getBlock(block)->getDataNonConst(j);
172 if ( !vecBCType_.at(loc).compare("Dirichlet") ) {
173 for (int dd=0; dd < dofs; dd++)
174 values[ dofs * i + dd ] = result.at(dd);
175 }
176 else if ( !vecBCType_.at(loc).compare("Dirichlet_X") ) {
177 values[ dofs * i + 0 ] = result.at(0);
178 }
179 else if ( !vecBCType_.at(loc).compare("Dirichlet_Y") ) {
180 values[ dofs * i + 1 ] = result.at(1);
181 }
182 else if ( !vecBCType_.at(loc).compare("Dirichlet_Z") ) {
183 values[ dofs * i + 2 ] = result.at(2);
184 }
185 else if ( !vecBCType_.at(loc).compare("Dirichlet_X_Y") ) {
186 values[ dofs * i + 0 ] = result.at(0);
187 values[ dofs * i + 1 ] = result.at(1);
188 }
189 else if ( !vecBCType_.at(loc).compare("Dirichlet_X_Z") ) {
190 values[ dofs * i + 0 ] = result.at(0);
191 values[ dofs * i + 2 ] = result.at(2);
193 else if ( !vecBCType_.at(loc).compare("Dirichlet_Y_Z") ) {
194 values[ dofs * i + 1 ] = result.at(1);
195 values[ dofs * i + 2 ] = result.at(2);
197 }
198 }
199 }
200 }
201 }
202 }
203 //Code below is experimental
204 for (int i=0; i<vecBCType_.size(); i++) {
205 if( vecBCType_.at(i) == "Neumann" && vecBlockID_.at(i)==block){
206 DomainPtr_Type domain = vecDomain_.at(i);
207 FEFacPtr_Type feFactory = Teuchos::rcp( new FEFac_Type() );
208 feFactory->addFE(domain);
209
210 std::vector<SC> funcParameter(3 , 0.);//0: order, 1:time, 2:surface flag (no specific flag is used in function above)
211 funcParameter[1] = t;
212 funcParameter[2] = vecFlag_.at(i);
213 int dim = domain->getDimension();
214 int dofs = vecDofs_.at(i);
215 MultiVectorPtr_Type a;
216 MultiVectorPtr_Type aUnique;
217 if (dofs>1){
218 a = Teuchos::rcp(new MultiVector_Type ( domain->getMapVecFieldRepeated() ) );
219 aUnique = Teuchos::rcp(new MultiVector_Type ( domain->getMapVecFieldUnique() ) );
220 feFactory->assemblySurfaceIntegralFlag( dim, domain->getFEType(), a, "Vector", vecBC_func_.at(i), funcParameter );
221 }
222 else{
223 a = Teuchos::rcp(new MultiVector_Type ( domain->getMapRepeated() ) );
224 aUnique = Teuchos::rcp(new MultiVector_Type ( domain->getMapUnique() ) );
225 feFactory->assemblySurfaceIntegralFlag( dim, domain->getFEType(), a, "Scalar", vecBC_func_.at(i), funcParameter );
226 }
227 aUnique->exportFromVector( a, false, "Add" );
228 MultiVectorPtr_Type mv = blockMV->getBlockNonConst(block);
229 mv->update(1.,*aUnique,1.);
230 }
231 }
232 }
233}
234
235template<class SC,class LO,class GO,class NO>
236void BCBuilder<SC,LO,GO,NO>::setDirichletBoundaryFromExternal(Teuchos::ArrayRCP<SC>& values/*values will be set to this vector*/, LO index, int loc, double time, std::string type, Teuchos::ArrayRCP<SC> valuesSubstract) const{
237
238 Teuchos::ArrayRCP<SC> valuesEx = vecExternalSol_[loc]->getDataNonConst(0);//We assume that our Multivector only has one column
239
240 (*pointPtr_)[0] = valuesEx[index]; // laplace solution in current point
241
242 vecBC_func_.at(loc)( &((*pointPtr_)[0]), &((*resultPtr_)[0]), time, &(vecBC_Parameters_.at(loc).at(0)));
243
244 int dofs = resultPtr_->size();
245 if ( !vecBCType_.at(loc).compare("Dirichlet") ) {
246 if (type == "standard") {
247 for (int dd=0; dd < dofs; dd++){
248 values[ dofs * index + dd ] = (*resultPtr_)[dd];
249 }
250 }
251 else if (type == "BCMinusVec"){
252 for (int dd=0; dd < dofs; dd++)
253 values[ dofs * index + dd ] = (*resultPtr_)[dd] - valuesSubstract[ dofs * index + dd ];
254 }
255 else if (type == "VecMinusBC"){
256 for (int dd=0; dd < dofs; dd++)
257 values[ dofs * index + dd ] = valuesSubstract[ dofs * index + dd ] - (*resultPtr_)[dd];
258 }
259 }
260 else {
261 TEUCHOS_TEST_FOR_EXCEPTION( true, std::runtime_error, "BCBuilder::setDirichletBoundaryFromExternal only implemented for full Dirichlet.");
262 }
263}
264
265template<class SC,class LO,class GO,class NO>
266void BCBuilder<SC,LO,GO,NO>::setBCMinusVector(const BlockMultiVectorPtr_Type &outBlockMV, const BlockMultiVectorPtr_Type &substractBlockMV, double t ) const{
267
268 TEUCHOS_TEST_FOR_EXCEPTION( outBlockMV->getNumVectors()>1, std::runtime_error, "BCBuilder::setRHS() only for getNumVectors == 1.");
269
270 vec_dbl_Type result;
271 vec_dbl_Type point;
272 int loc = 0;
273
274 for (int block = 0; block < outBlockMV->size(); block++) { // blocks of RHS vector
275 if(blockHasDirichletBC(block, loc)){
276 if (vecFlowRateBool_[loc]) { // we use the external vector here with flowrate
277 this->determineVelocityForFlowrate(loc,t);
278 }
279 vec_int_ptr_Type bcFlags = vecDomain_.at(loc)->getBCFlagUnique();
280 vec2D_dbl_ptr_Type vecPoints = vecDomain_.at(loc)->getPointsUnique();
281 int dim = vecDomain_.at(loc)->getDimension();
282 int dofs = vecDofs_.at(loc);
283
284 result.resize(dofs,0.);
285 point.resize(dim,0.);
286
287 resultPtr_ = Teuchos::rcp(new vec_dbl_Type(dofs,0.));
288 pointPtr_ = Teuchos::rcp(new vec_dbl_Type(dim,0.));
289
290 Teuchos::ArrayRCP<SC> valuesRHS = outBlockMV->getBlock(block)->getDataNonConst(0);
291 Teuchos::ArrayRCP<SC> valuesSubstract = substractBlockMV->getBlock(block)->getDataNonConst(0);
292
293 for (int i = 0 ; i < bcFlags->size(); i++) {
294 int flag = bcFlags->at(i);
295 if(findFlag(flag,block, loc)){ // loc changed to postion corresponding to bcFlag
296 if (!vecBCType_.at(loc).compare("Dirichlet") ||
297 !vecBCType_.at(loc).compare("Dirichlet_X") ||
298 !vecBCType_.at(loc).compare("Dirichlet_Y") ||
299 !vecBCType_.at(loc).compare("Dirichlet_Z")||
300 !vecBCType_.at(loc).compare("Dirichlet_X_Y") ||
301 !vecBCType_.at(loc).compare("Dirichlet_X_Z") ||
302 !vecBCType_.at(loc).compare("Dirichlet_Y_Z")) {
303 if (!vecExternalSol_[loc].is_null()) { // we use the external vector here
304 std::string type = "BCMinusVec";
305 this->setDirichletBoundaryFromExternal(valuesRHS, i, loc, t, type, valuesSubstract);
306 }
307 else{
308 for (int d=0; d < dim; d++) {
309 point.at(d) = vecPoints->at(i).at(d);
310 }
311 for (int d=0; d < dofs; d++) {
312 result.at(d) = vecPoints->at(i).at(d);
313 }
314 vecBC_func_.at(loc)(&point.at(0), &result.at(0), t, &(vecBC_Parameters_.at(loc).at(0)));
315
316 for (UN j=0; j<outBlockMV->getNumVectors(); j++) {
317 Teuchos::ArrayRCP<SC> values = outBlockMV->getBlock(block)->getDataNonConst(j);
318 Teuchos::ArrayRCP<const SC> substractValues = substractBlockMV->getBlock(block)->getData(j);
319 if ( !vecBCType_.at(loc).compare("Dirichlet") ) {
320 for (int dd=0; dd < dofs; dd++)
321 values[ dofs * i + dd ] = result.at(dd) - substractValues[ dofs * i + dd ];
322 }
323 else if ( !vecBCType_.at(loc).compare("Dirichlet_X") ) {
324 values[ dofs * i + 0 ] = result.at(0) - substractValues[ dofs * i + 0 ];
325 }
326 else if ( !vecBCType_.at(loc).compare("Dirichlet_Y") ) {
327 values[ dofs * i + 1 ] = result.at(1) - substractValues[ dofs * i + 1 ];
328 }
329 else if ( !vecBCType_.at(loc).compare("Dirichlet_Z") ) {
330 values[ dofs * i + 2 ] = result.at(2) - substractValues[ dofs * i + 2 ];
331 }
332 else if ( !vecBCType_.at(loc).compare("Dirichlet_X_Y") ) {
333 values[ dofs * i + 0 ] = result.at(0) - substractValues[ dofs * i + 0 ];
334 values[ dofs * i + 1 ] = result.at(1) - substractValues[ dofs * i + 1 ];
335 }
336 else if ( !vecBCType_.at(loc).compare("Dirichlet_X_Z") ) {
337 values[ dofs * i + 0 ] = result.at(0) - substractValues[ dofs * i + 0 ];
338 values[ dofs * i + 2 ] = result.at(2) - substractValues[ dofs * i + 2 ];
339 }
340 else if ( !vecBCType_.at(loc).compare("Dirichlet_Y_Z") ) {
341 values[ dofs * i + 1 ] = result.at(1) - substractValues[ dofs * i + 1 ];
342 values[ dofs * i + 2 ] = result.at(2) - substractValues[ dofs * i + 2 ];
343 }
344 }
345 }
346 }
347 }
348 }
349 }
350 //Code below is experimental
351 for (int i=0; i<vecBCType_.size(); i++) {
352 if( vecBCType_.at(i) == "Neumann" && vecBlockID_.at(i)==block){
353 DomainPtr_Type domain = vecDomain_.at(i);
354 FEFacPtr_Type feFactory = Teuchos::rcp( new FEFac_Type() );
355 feFactory->addFE(domain);
356
357 std::vector<SC> funcParameter(3 , 0.);//0: order, 1:time, 2:surface flag (no specific flag is used in function above)
358 funcParameter[1] = t;
359 funcParameter[2] = vecFlag_.at(i);
360 int dim = domain->getDimension();
361 int dofs = vecDofs_.at(i);
362 MultiVectorPtr_Type a;
363 MultiVectorPtr_Type aUnique;
364 if (dofs>1){
365 a = Teuchos::rcp(new MultiVector_Type ( domain->getMapVecFieldRepeated() ) );
366 aUnique = Teuchos::rcp(new MultiVector_Type ( domain->getMapVecFieldUnique() ) );
367 feFactory->assemblySurfaceIntegralFlag( dim, domain->getFEType(), a, "Vector", vecBC_func_.at(i), funcParameter );
368 }
369 else{
370 a = Teuchos::rcp(new MultiVector_Type ( domain->getMapRepeated() ) );
371 aUnique = Teuchos::rcp(new MultiVector_Type ( domain->getMapUnique() ) );
372 feFactory->assemblySurfaceIntegralFlag( dim, domain->getFEType(), a, "Scalar", vecBC_func_.at(i), funcParameter );
373 }
374 aUnique->exportFromVector( a, false, "Add" );
375 MultiVectorPtr_Type mv = outBlockMV->getBlockNonConst(block);
376 mv->update(1.,*aUnique,1.);
377 MultiVectorPtr_Type mvMinus = substractBlockMV->getBlockNonConst(block);
378 mv->update(-1.,*mvMinus,1.);
379 }
380 }
381 }
382}
383
384template<class SC,class LO,class GO,class NO>
385void BCBuilder<SC,LO,GO,NO>::setVectorMinusBC(const BlockMultiVectorPtr_Type &outBlockMV, const BlockMultiVectorPtr_Type &substractBlockMV, double t) const{
386 vec_dbl_Type result;
387 vec_dbl_Type point;
388 int loc = 0;
389
390 TEUCHOS_TEST_FOR_EXCEPTION( outBlockMV->getNumVectors()>1, std::runtime_error, "BCBuilder::setRHS() only for getNumVectors == 1.");
391
392 for (int block = 0; block < outBlockMV->size(); block++) { // blocks of RHS vector
393 if(blockHasDirichletBC(block, loc)){
394 if (vecFlowRateBool_[loc]) { // we use the external vector here with flowrate
395 this->determineVelocityForFlowrate(loc,t);
396 }
397 vec_int_ptr_Type bcFlags = vecDomain_.at(loc)->getBCFlagUnique();
398 vec2D_dbl_ptr_Type vecPoints = vecDomain_.at(loc)->getPointsUnique();
399 int dim = vecDomain_.at(loc)->getDimension();
400 int dofs = vecDofs_.at(loc);
401
402 result.resize(dofs,0.);
403 point.resize(dim,0.);
404
405 resultPtr_ = Teuchos::rcp(new vec_dbl_Type(dofs,0.));
406 pointPtr_ = Teuchos::rcp(new vec_dbl_Type(dim,0.));
407
408 Teuchos::ArrayRCP<SC> valuesRHS = outBlockMV->getBlock(block)->getDataNonConst(0);
409 Teuchos::ArrayRCP<SC> valuesSubstract = substractBlockMV->getBlock(block)->getDataNonConst(0);
410 for (int i = 0 ; i < bcFlags->size(); i++) {
411 int flag = bcFlags->at(i);
412 if(findFlag(flag,block, loc)){ // loc changed to postion corresponding to bcFlag
413 if (!vecBCType_.at(loc).compare("Dirichlet") ||
414 !vecBCType_.at(loc).compare("Dirichlet_X") ||
415 !vecBCType_.at(loc).compare("Dirichlet_Y") ||
416 !vecBCType_.at(loc).compare("Dirichlet_Z") ||
417 !vecBCType_.at(loc).compare("Dirichlet_X_Y") ||
418 !vecBCType_.at(loc).compare("Dirichlet_X_Z") ||
419 !vecBCType_.at(loc).compare("Dirichlet_Y_Z")) {
420 if (!vecExternalSol_[loc].is_null()) { // we use the external vector here
421 std::string type = "VecMinusBC";
422 this->setDirichletBoundaryFromExternal(valuesRHS, i, loc, t, type, valuesSubstract);
423 }
424 else{
425 for (int d=0; d < dim; d++) {
426 point.at(d) = vecPoints->at(i).at(d);
427 }
428 for (int d=0; d < dofs; d++) {
429 result.at(d) = vecPoints->at(i).at(d);
430 }
431 vecBC_func_.at(loc)(&point.at(0), &result.at(0), t, &(vecBC_Parameters_.at(loc).at(0)));
432
433 for (UN j=0; j<outBlockMV->getNumVectors(); j++) {
434 Teuchos::ArrayRCP<SC> values = outBlockMV->getBlock(block)->getDataNonConst(j);
435 Teuchos::ArrayRCP<const SC> substractValues = substractBlockMV->getBlock(block)->getData(j);
436 if ( !vecBCType_.at(loc).compare("Dirichlet") ) {
437 for (int dd=0; dd < dofs; dd++)
438 values[ dofs * i + dd ] = substractValues[ dofs * i + dd ] - result.at(dd);
439 }
440 else if ( !vecBCType_.at(loc).compare("Dirichlet_X") ) {
441 values[ dofs * i + 0 ] = substractValues[ dofs * i + 0 ] - result.at(0);
442 }
443 else if ( !vecBCType_.at(loc).compare("Dirichlet_Y") ) {
444 values[ dofs * i + 1 ] = substractValues[ dofs * i + 1 ] - result.at(1);
445 }
446 else if ( !vecBCType_.at(loc).compare("Dirichlet_Z") ) {
447 values[ dofs * i + 2 ] = substractValues[ dofs * i + 2 ] - result.at(2);
448 }
449 else if ( !vecBCType_.at(loc).compare("Dirichlet_X_Y") ) {
450 values[ dofs * i + 0 ] = substractValues[ dofs * i + 0 ] - result.at(0);
451 values[ dofs * i + 1 ] = substractValues[ dofs * i + 1 ] - result.at(1);
452 }
453 else if ( !vecBCType_.at(loc).compare("Dirichlet_X_Z") ) {
454 values[ dofs * i + 0 ] = substractValues[ dofs * i + 0 ] - result.at(0);
455 values[ dofs * i + 2 ] = substractValues[ dofs * i + 2 ] - result.at(2);
456 }
457 else if ( !vecBCType_.at(loc).compare("Dirichlet_Y_Z") ) {
458 values[ dofs * i + 1 ] = substractValues[ dofs * i + 1 ] - result.at(1);
459 values[ dofs * i + 2 ] = substractValues[ dofs * i + 2 ] - result.at(2);
460 }
461 }
462 }
463 }
464 }
465 }
466 }
467 //Code below is experimental
468 for (int i=0; i<vecBCType_.size(); i++) {
469 if( vecBCType_.at(i) == "Neumann" && vecBlockID_.at(i)==block){
470 DomainPtr_Type domain = vecDomain_.at(i);
471 FEFacPtr_Type feFactory = Teuchos::rcp( new FEFac_Type() );
472 feFactory->addFE(domain);
473
474 std::vector<SC> funcParameter(3 , 0.);//0: order, 1:time, 2:surface flag (no specific flag is used in function above)
475 funcParameter[1] = t;
476 funcParameter[2] = vecFlag_.at(i);
477 int dim = domain->getDimension();
478 int dofs = vecDofs_.at(i);
479 MultiVectorPtr_Type a;
480 MultiVectorPtr_Type aUnique;
481 if (dofs>1){
482 a = Teuchos::rcp(new MultiVector_Type ( domain->getMapVecFieldRepeated() ) );
483 aUnique = Teuchos::rcp(new MultiVector_Type ( domain->getMapVecFieldUnique() ) );
484 feFactory->assemblySurfaceIntegralFlag( dim, domain->getFEType(), a, "Vector", vecBC_func_.at(i), funcParameter );
485 }
486 else{
487 a = Teuchos::rcp(new MultiVector_Type ( domain->getMapRepeated() ) );
488 aUnique = Teuchos::rcp(new MultiVector_Type ( domain->getMapUnique() ) );
489 feFactory->assemblySurfaceIntegralFlag( dim, domain->getFEType(), a, "Scalar", vecBC_func_.at(i), funcParameter );
490 }
491 aUnique->exportFromVector( a, false, "Add" );
492 MultiVectorPtr_Type mv = outBlockMV->getBlockNonConst(block);
493 mv->update(1.,*aUnique,1.);
494 MultiVectorPtr_Type mvMinus = substractBlockMV->getBlockNonConst(block);
495 mv->update(1.,*mvMinus,-1.);
496 }
497 }
498 }
499}
500
501template<class SC,class LO,class GO,class NO>
502void BCBuilder<SC,LO,GO,NO>::setAllDirichletZero(const BlockMultiVectorPtr_Type &blockMV) const{
503 vec_dbl_Type result;
504 int loc = 0;
505 for (int block = 0; block < blockMV->size(); block++) { // blocks of RHS vector
506 if(blockHasDirichletBC(block, loc)){
507 vec_int_ptr_Type bcFlags = vecDomain_.at(loc)->getBCFlagUnique();
508 vec2D_dbl_ptr_Type vecPoints = vecDomain_.at(loc)->getPointsUnique();
509 int dim = vecDomain_.at(loc)->getDimension();
510 int dofs = vecDofs_.at(loc);
511
512 result.resize(dofs,0.);
513
514 for (int i = 0 ; i < bcFlags->size(); i++) {
515 int flag = bcFlags->at(i);
516 if(findFlag(flag,block, loc)){ // loc changed to postion corresponding to bcFlag
517 if ( !vecBCType_.at(loc).compare("Dirichlet") ||
518 !vecBCType_.at(loc).compare("Dirichlet_X") ||
519 !vecBCType_.at(loc).compare("Dirichlet_Y") ||
520 !vecBCType_.at(loc).compare("Dirichlet_Z") ||
521 !vecBCType_.at(loc).compare("Dirichlet_X_Y") ||
522 !vecBCType_.at(loc).compare("Dirichlet_X_Z") ||
523 !vecBCType_.at(loc).compare("Dirichlet_Y_Z")
524 ){
525
526 for (UN j=0; j<blockMV->getNumVectors(); j++) {
527 Teuchos::ArrayRCP<SC> values = blockMV->getBlock(block)->getDataNonConst(j);
528 if ( !vecBCType_.at(loc).compare("Dirichlet") ) {
529 for (int dd=0; dd < dofs; dd++)
530 values[ dofs * i + dd ] = result[dd];
531 }
532 else if ( !vecBCType_.at(loc).compare("Dirichlet_X") ) {
533 values[ dofs * i + 0 ] = result[0];
534 }
535 else if ( !vecBCType_.at(loc).compare("Dirichlet_Y") ) {
536 values[ dofs * i + 1 ] = result[1];
537 }
538 else if ( !vecBCType_.at(loc).compare("Dirichlet_Z") ) {
539 values[ dofs * i + 2 ] = result[2];
540 }
541 else if ( !vecBCType_.at(loc).compare("Dirichlet_X_Y") ) {
542 values[ dofs * i + 0 ] = result[0];
543 values[ dofs * i + 1 ] = result[1];
544 }
545 else if ( !vecBCType_.at(loc).compare("Dirichlet_X_Z") ) {
546 values[ dofs * i + 0 ] = result[0];
547 values[ dofs * i + 2 ] = result[2];
548 }
549 else if ( !vecBCType_.at(loc).compare("Dirichlet_Y_Z") ) {
550 values[ dofs * i + 1 ] = result[1];
551 values[ dofs * i + 2 ] = result[2];
552 }
553 }
554 }
555 }
556 }
557 }
558 }
559}
560
561// We have the problem, that we want to prescribe a constant flow rate
562// As a boundary condition, we need to prescibe a certain velocity.
563// In an FSI setting, the flowrate is also influences by the changing
564// area of the inlet. Consequently, the velocity we prescibe changes depending
565// on the changing area.
566// We follow the following idea: The flowrate Q is given as \int_{Inlet} u * n dA
567// In our case we have a parabolic-like inflow profile, which we can write as u = u_max * vec
568// With the desired flow rate Q, we can determine u_max as follows:
569// \int_{Inlet} u_max * vec * n dA != Q
570// <=> u_max * \int_{Inlet} vec * n dA != Q
571// <=> u_max = Q / \int_{Inlet} vec * n dA
572template<class SC,class LO,class GO,class NO>
573void BCBuilder<SC,LO,GO,NO>::determineVelocityForFlowrate(LO i, double time) const{
574
575 // Domain of the corresponing boundary condition. Most probably fluid domain
576 DomainPtr_Type domain = vecDomain_.at(i);
577 // A vector containing and desribing the parabolic inflow
578 MultiVectorConstPtr_Type parabolic_unique_const = vecExternalSol_[i]; // e.g. Laplace soution on inlet
579 MultiVectorPtr_Type parabolic_unique = Teuchos::rcp_const_cast<MultiVector_Type> ( parabolic_unique_const ); // dirty const casting
580
581 MultiVectorPtr_Type parabolic_rep = Teuchos::rcp(new MultiVector_Type ( domain->getMapRepeated() ) );
582 // We normalize the solution and distribute it to the repeated map
583 SC maxValue = parabolic_unique->getMax();
584 parabolic_unique->scale(1./maxValue); // normalizing solution
585 parabolic_rep->importFromVector(parabolic_unique,false,"Insert");
586 // We define an FE Factory
587 FEFacPtr_Type feFactory = Teuchos::rcp( new FEFac_Type() );
588 feFactory->addFE(domain);
589 // We determine the current desired flow rate via the input function
590 std::vector<SC> funcParameter = vecBC_Parameters_[i];
591 vec_dbl_Type p1 = {1.,1.,1.}; // Dummy vector
592 vec_dbl_Type flowRate = {0.}; //
593 vecBC_func_flowRate_.at(i)( &(p1[0]), &(flowRate[0]), time, &(funcParameter[0])); // Determine Flowrate based on BC function
594 // We assemble the flowrate for parabolic inflow profile. As we have the inflow profile normalized we have: vec* u_max = RB Inlet normally
595 double flowRateParabolic=0.;
596 feFactory->assemblyFlowRate(domain->getDimension(), flowRateParabolic, domain->getFEType(),1, vecFlag_[i] , parabolic_rep);
597 // Then we have \int_{Inlet} vec * n dA * u_max == Q <=> u_max = Q/\int_{Inlet} vec * n dA, and Q is given as 'desired flowrate' in flowrate
598 double maxVelocity = flowRate[0] / std::fabs(flowRateParabolic);
599 // Then we replace the parameter which contains the maximum velocity with the updated one
600 vecBC_Parameters_[i][0] = maxVelocity;
601
602}
603
604
605
606template<class SC,class LO,class GO,class NO>
608 int dummyloc;
609
610 return blockHasDirichletBC(block, dummyloc);
611}
612
613template<class SC,class LO,class GO,class NO>
614bool BCBuilder<SC,LO,GO,NO>::blockHasDirichletBC(int block, int &loc) const{
615
616 bool hasBC = false;
617 bool found = false;
618 loc = -1;
619 std::vector<int>::const_iterator it = vecBlockID_.begin();
620 while (it!=vecBlockID_.end() && !found) {
621 it = std::find(it,vecBlockID_.end(),block);
622 if (it!=vecBlockID_.end()) {
623 loc = distance(vecBlockID_.begin(),it);
624 it++;
625 if (!vecBCType_.at(loc).compare("Dirichlet") ||
626 !vecBCType_.at(loc).compare("Dirichlet_X") ||
627 !vecBCType_.at(loc).compare("Dirichlet_Y") ||
628 !vecBCType_.at(loc).compare("Dirichlet_Z") ||
629 !vecBCType_.at(loc).compare("Dirichlet_X_Y") ||
630 !vecBCType_.at(loc).compare("Dirichlet_X_Z") ||
631 !vecBCType_.at(loc).compare("Dirichlet_Y_Z")
632 ) {
633 found = true;
634 hasBC = true;
635 }
636 }
637 }
638
639 return hasBC;
640}
641
642template<class SC,class LO,class GO,class NO>
643bool BCBuilder<SC,LO,GO,NO>::findFlag(LO flag, int block, int &loc) const{
644
645#ifdef BCBuilder_TIMER
646 Teuchos::TimeMonitor FindFlagMonitor(*FindFlagTimer_);
647#endif
648
649 bool hasFlag = false;
650 bool found = false;
651 loc = -1;
652 std::vector<int>::const_iterator it = vecFlag_.begin();
653 while (it!=vecFlag_.end() && !found) {
654 it = std::find(it,vecFlag_.end(),flag);
655 if (it!=vecFlag_.end()) {
656 loc = distance(vecFlag_.begin(),it);
657 if (vecBlockID_.at(loc)==block) {
658 hasFlag = true;
659 found = true;
660 }
661 else{
662 it++;
663 }
664 }
665 }
666
667 return hasFlag;
668}
669
670
671template<class SC,class LO,class GO,class NO>
672void BCBuilder<SC,LO,GO,NO>::setSystemScaled(const BlockMatrixPtr_Type &blockMatrix,double eps) const{
673
674 UN numBlocks = blockMatrix->size();
675 int loc;
676 for (UN blockRow = 0; blockRow < numBlocks; blockRow++) {
677
678#ifdef BCBuilder_TIMER
679 Teuchos::TimeMonitor SetSystemRowMonitor(*SetSystemRowTimer_);
680#endif
681 bool boolBlockHasDirichlet = false;
682
683 {
684#ifdef BCBuilder_TIMER
685 Teuchos::TimeMonitor BlockRowHasDirichletMonitor(*BlockRowHasDirichletTimer_);
686#endif
687 boolBlockHasDirichlet = blockHasDirichletBC(blockRow,loc);
688 }
689
690 if (boolBlockHasDirichlet){
691 for (UN blockCol = 0; blockCol < numBlocks ; blockCol++) {
692 if ( blockMatrix->blockExists( blockRow, blockCol ) ) {
693 MatrixPtr_Type matrix = blockMatrix->getBlock( blockRow, blockCol );
694 setDirichletBCScaled( matrix, loc, blockRow, blockRow==blockCol,eps );
695 }
696 }
697 }
698 }
699}
700
701
702template<class SC,class LO,class GO,class NO>
703void BCBuilder<SC,LO,GO,NO>::setSystem(const BlockMatrixPtr_Type &blockMatrix) const{
704
705 UN numBlocks = blockMatrix->size();
706 int loc;
707 for (UN blockRow = 0; blockRow < numBlocks; blockRow++) {
708
709#ifdef BCBuilder_TIMER
710 Teuchos::TimeMonitor SetSystemRowMonitor(*SetSystemRowTimer_);
711#endif
712 bool boolBlockHasDirichlet = false;
713
714 {
715#ifdef BCBuilder_TIMER
716 Teuchos::TimeMonitor BlockRowHasDirichletMonitor(*BlockRowHasDirichletTimer_);
717#endif
718
719 boolBlockHasDirichlet = blockHasDirichletBC(blockRow,loc);
720 }
721
722 if (boolBlockHasDirichlet){
723 for (UN blockCol = 0; blockCol < numBlocks ; blockCol++) {
724 if (blockMatrix->blockExists(blockRow, blockCol)) {
725 MatrixPtr_Type matrix = blockMatrix->getBlock(blockRow, blockCol);
726 setDirichletBC(matrix, loc, blockRow, blockRow == blockCol);
727 } else if (blockRow == blockCol) {
728 // If the current block is a diagonal block but does not exist we need to create it and insert the
729 // value 1 in all diagonal entries of the block corresponding to a Dirichlet boundary condition. We
730 // only need max. one entry per row.
731
732 // To get the correct domain from vecDomain_ we need to know the index of any (in this case the
733 // first) entry that was done with addBC for this block
734 auto matrixMap = this->vecDomain_.at(loc)->getMapUnique();
735 // // Use Xpetra::MatrixFactory to build a matrix with known row and column maps
736 auto tpetraMatrix = Teuchos::rcp(new Tpetra::CrsMatrix<SC,LO,GO,NO>(matrixMap->getTpetraMap(), matrixMap->getTpetraMap(), 1)); // Tpetra::MatrixFactory<SC, LO, GO, NO>::Build(matrixMap->getTpetraMap(), matrixMap->getTpetraMap(), 1);
737
738 MatrixPtr_Type matrix = Teuchos::rcp(new Matrix_Type(tpetraMatrix));
739
740 // Fill the matrix with zeros on the diagonal
741 Teuchos::Array<LO> colIndex(1, 0);
742 Teuchos::Array<SC> val(1, Teuchos::ScalarTraits<SC>::zero());
743 for (auto i = 0; i < matrixMap->getNodeNumElements(); i++) {
744 colIndex[0] = i;
745 // Since the matrix has a column map we can use insertLocalValues(). This is more efficient!
746 matrix->insertLocalValues(i, colIndex, val);
747 }
748 matrix->fillComplete(matrixMap, matrixMap);
749 setDirichletBC(matrix, loc, blockRow, true);
750 blockMatrix->addBlock(matrix, blockRow, blockCol);
751 }
752 }
753 }
754 }
755}
756
757
758template<class SC,class LO,class GO,class NO>
759void BCBuilder<SC,LO,GO,NO>::setDirichletBCScaled(const MatrixPtr_Type &matrix, int loc, int blockRow, bool isDiagonalBlock, double eps) const{
760
761 matrix->resumeFill();
762 bool isDirichlet;
763 UN dofsPerNode = vecDofs_.at(loc);
764 vec_int_ptr_Type bcFlags = vecDomain_.at(loc)->getBCFlagUnique();
765 MapConstPtr_Type nodeMap = vecDomain_.at(loc)->getMapUnique();
766
767 for (LO i=0; i<bcFlags->size(); i++) { // flags are for nodes.
768 int flag = bcFlags->at(i);
769 isDirichlet = false;
770 int locThisFlag;
771
772 if( findFlag(flag, blockRow, locThisFlag) ){
773
774 if( !vecBCType_.at(locThisFlag).compare("Dirichlet") ||
775 !vecBCType_.at(locThisFlag).compare("Dirichlet_X") ||
776 !vecBCType_.at(locThisFlag).compare("Dirichlet_Y") ||
777 !vecBCType_.at(locThisFlag).compare("Dirichlet_Z") ||
778 !vecBCType_.at(locThisFlag).compare("Dirichlet_X_Y") ||
779 !vecBCType_.at(locThisFlag).compare("Dirichlet_X_Z") ||
780 !vecBCType_.at(locThisFlag).compare("Dirichlet_Y_Z") )
781 isDirichlet = true;
782 if (isDirichlet) {
783
784 if (isDiagonalBlock)
785 setLocalRowEntry(matrix, i, dofsPerNode, locThisFlag,eps );
786 else
787 setLocalRowZero(matrix, i, dofsPerNode, locThisFlag );
788 }
789 }
790 }
791 matrix->fillComplete( matrix->getMap("domain"), matrix->getMap("range") );
792}
793
794
795template<class SC,class LO,class GO,class NO>
796void BCBuilder<SC,LO,GO,NO>::setLocalRowEntry(const MatrixPtr_Type &matrix, LO localNode, UN dofsPerNode, int loc,double eps) const{
797
798 Teuchos::ArrayView<const SC> valuesOld;
799 Teuchos::ArrayView<const LO> indices;
800
801 MapConstPtr_Type colMap = matrix->getMap("col");
802 LO localDof = (LO) (dofsPerNode * localNode);
803 for (UN dof=0; dof<dofsPerNode; dof++) {
804 if ( vecBCType_.at(loc) == "Dirichlet" ||
805 (vecBCType_.at(loc) == "Dirichlet_X" && dof==0 ) ||
806 (vecBCType_.at(loc) == "Dirichlet_Y" && dof==1 ) ||
807 (vecBCType_.at(loc) == "Dirichlet_Z" && dof==2 ) ||
808 (vecBCType_.at(loc) == "Dirichlet_X_Y" && dof!=2 )||
809 (vecBCType_.at(loc) == "Dirichlet_X_Z" && dof!=1 )||
810 (vecBCType_.at(loc) == "Dirichlet_Y_Z" && dof!=0 )
811 ) {
812 // cout << " Setting Dirichlet Row 1 for node " << localDof << " of type " << vecBCType_.at(loc) <<endl;
813 GO globalDof = matrix->getMap()->getGlobalElement( localDof );
814 matrix->getLocalRowView(localDof, indices, valuesOld);
815 double rowSum = 0.;
816 for (UN j=0; j<indices.size(); j++) {
817 rowSum += abs(valuesOld[j]);
818 }
819 Teuchos::Array<SC> values( valuesOld.size(), Teuchos::ScalarTraits<SC>::zero() );
820 bool setOne = false;
821 for (UN j=0; j<indices.size() && !setOne; j++) {
822 if ( colMap->getGlobalElement( indices[j] ) == globalDof ){
823 values[j] = valuesOld[j]*eps;
824
825 setOne = true;
826 }
827 }
828 matrix->replaceLocalValues(localDof, indices(), values());
829 }
830 localDof++;
831 }
832
833}
834
835
836template<class SC,class LO,class GO,class NO>
837void BCBuilder<SC,LO,GO,NO>::setDirichletBC(const MatrixPtr_Type &matrix, int loc, int blockRow, bool isDiagonalBlock) const{
838
839 matrix->resumeFill();
840 bool isDirichlet;
841 UN dofsPerNode = vecDofs_.at(loc);
842 vec_int_ptr_Type bcFlags = vecDomain_.at(loc)->getBCFlagUnique();
843
844 for (LO i=0; i<bcFlags->size(); i++) { // flags are for nodes.
845 int flag = bcFlags->at(i);
846 isDirichlet = false;
847 int locThisFlag;
848
849 if( findFlag(flag, blockRow, locThisFlag) ){
850
851 if( !vecBCType_.at(locThisFlag).compare("Dirichlet") ||
852 !vecBCType_.at(locThisFlag).compare("Dirichlet_X") ||
853 !vecBCType_.at(locThisFlag).compare("Dirichlet_Y") ||
854 !vecBCType_.at(locThisFlag).compare("Dirichlet_Z") ||
855 !vecBCType_.at(locThisFlag).compare("Dirichlet_X_Y") ||
856 !vecBCType_.at(locThisFlag).compare("Dirichlet_X_Z") ||
857 !vecBCType_.at(locThisFlag).compare("Dirichlet_Y_Z") )
858 isDirichlet = true;
859 if (isDirichlet) {
860
861 if (isDiagonalBlock)
862 setLocalRowOne(matrix, i, dofsPerNode, locThisFlag );
863 else
864 setLocalRowZero(matrix, i, dofsPerNode, locThisFlag );
865 }
866 }
867 }
868 matrix->fillComplete( matrix->getMap("domain"), matrix->getMap("range") );
869}
870
871template<class SC,class LO,class GO,class NO>
872void BCBuilder<SC,LO,GO,NO>::setLocalRowOne(const MatrixPtr_Type &matrix, LO localNode, UN dofsPerNode, int loc) const{
873
874 Teuchos::ArrayView<const SC> valuesOld;
875 Teuchos::ArrayView<const LO> indices;
876
877 MapConstPtr_Type colMap = matrix->getMap("col");
878 LO localDof = (LO) (dofsPerNode * localNode);
879 for (UN dof=0; dof<dofsPerNode; dof++) {
880 if ( vecBCType_.at(loc) == "Dirichlet" ||
881 (vecBCType_.at(loc) == "Dirichlet_X" && dof==0 ) ||
882 (vecBCType_.at(loc) == "Dirichlet_Y" && dof==1 ) ||
883 (vecBCType_.at(loc) == "Dirichlet_Z" && dof==2 ) ||
884 (vecBCType_.at(loc) == "Dirichlet_X_Y" && dof!=2 )||
885 (vecBCType_.at(loc) == "Dirichlet_X_Z" && dof!=1 )||
886 (vecBCType_.at(loc) == "Dirichlet_Y_Z" && dof!=0 )
887 ) {
888 // std::cout << " Setting Dirichlet Row 1 for node " << localDof << " of type " << vecBCType_.at(loc) <<endl;
889 GO globalDof = matrix->getMap()->getGlobalElement( localDof );
890 matrix->getLocalRowView(localDof, indices, valuesOld);
891 Teuchos::Array<SC> values( valuesOld.size(), Teuchos::ScalarTraits<SC>::zero() );
892 bool setOne = false;
893 for (UN j=0; j<indices.size() && !setOne; j++) {
894 if ( colMap->getGlobalElement( indices[j] ) == globalDof ){
895 values[j] = Teuchos::ScalarTraits<SC>::one();
896 setOne = true;
897 }
898 }
899 matrix->replaceLocalValues(localDof, indices(), values());
900 }
901 localDof++;
902 }
903
904}
905
906template<class SC,class LO,class GO,class NO>
907void BCBuilder<SC,LO,GO,NO>::setLocalRowZero(const MatrixPtr_Type &matrix, LO localNode, UN dofsPerNode, int loc) const{
908
909 Teuchos::ArrayView<const SC> valuesOld;
910 Teuchos::ArrayView<const LO> indices;
911 LO localDof = (LO) (dofsPerNode * localNode);
912 for (UN dof=0; dof<dofsPerNode; dof++) {
913 if ( vecBCType_.at(loc) == "Dirichlet" ||
914 (vecBCType_.at(loc) == "Dirichlet_X" && dof==0 ) ||
915 (vecBCType_.at(loc) == "Dirichlet_Y" && dof==1 ) ||
916 (vecBCType_.at(loc) == "Dirichlet_Z" && dof==2 ) ||
917 (vecBCType_.at(loc) == "Dirichlet_X_Y" && dof!=2 )||
918 (vecBCType_.at(loc) == "Dirichlet_X_Z" && dof!=1 )||
919 (vecBCType_.at(loc) == "Dirichlet_Y_Z" && dof!=0 )
920 ) {
921 matrix->getLocalRowView(localDof, indices, valuesOld);
922 Teuchos::Array<SC> values( valuesOld.size(), Teuchos::ScalarTraits<SC>::zero() );
923 matrix->replaceLocalValues(localDof, indices(), values());
924 }
925 localDof++;
926 }
927}
928
929
930
931template<class SC,class LO,class GO,class NO>
933 std::vector<int>::iterator it;
934 it = find(vecBlockID_.begin(),vecBlockID_.end(),block);
935#ifdef ASSERTS_WARNINGS
936 MYASSERT(it!=vecBlockID_.end(),"Requested dof information is not known to BCBuilder class.");
937#endif
938
939 return vecDofs_.at(distance(vecBlockID_.begin(),it));
940}
941}
942
943#endif
void setBCMinusVector(const BlockMultiVectorPtr_Type &outBlockMV, const BlockMultiVectorPtr_Type &substractBlockMV, double t=0.) const
Definition BCBuilder_def.hpp:266
void setVectorMinusBC(const BlockMultiVectorPtr_Type &outBlockMV, const BlockMultiVectorPtr_Type &substractBlockMV, double t=0.) const
Definition BCBuilder_def.hpp:385
void setRHS(const BlockMultiVectorPtr_Type &blockMV, double t=0.) const
void set(const BlockMatrixPtr_Type &blockMatrix, const BlockMultiVectorPtr_Type &b, double t=0.) const
Setting bundary condtions to problem.
Definition BCBuilder_def.hpp:114
void setSystem(const BlockMatrixPtr_Type &blockMatrix) const
void addBC(BC_func_Type funcBC, int flag, int block, const DomainPtr_Type &domain, std::string type, int dofs)
Adding Boundary Condition.
Definition BCBuilder_def.hpp:46
void setLocalRowZero(const MatrixPtr_Type &matrix, LO localNode, UN dofsPerNode, int loc) const
Definition BCBuilder_def.hpp:907
void setDirichletBoundaryFromExternal(Teuchos::ArrayRCP< SC > &values, LO index, int loc, double time, std::string type, Teuchos::ArrayRCP< SC > valuesSubstract=Teuchos::null) const
void setSystemScaled(const BlockMatrixPtr_Type &blockMatrix, double eps=1.0) const
Set boundary conditions to system.
Definition BCBuilder_def.hpp:672
void setDirichletBC(const MatrixPtr_Type &matrix, int loc, int blockRow, bool isDiagonalBlock) const
Definition BCBuilder_def.hpp:837
void setAllDirichletZero(const BlockMultiVectorPtr_Type &blockMV) const
Definition BCBuilder_def.hpp:502
bool blockHasDirichletBC(int block) const
void setDirichletBCScaled(const MatrixPtr_Type &matrix, int loc, int blockRow, bool isDiagonalBlock, double eps=1.0) const
Definition BCBuilder_def.hpp:759
bool findFlag(LO flag, int block, int &loc) const
void setLocalRowEntry(const MatrixPtr_Type &matrix, LO localNode, UN dofsPerNode, int loc, double eps) const
Definition BCBuilder_def.hpp:796
void setLocalRowOne(const MatrixPtr_Type &matrix, LO localNode, UN dofsPerNode, int loc) const
Definition BCBuilder_def.hpp:872
int dofsPerNodeAtBlock(int block)
Definition BCBuilder_def.hpp:932
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement_decl.hpp:36