00001
00002 #include <stdlib.h>
00003 #include <algorithm>
00004 #include <iostream>
00005
00006 #include "GLpApp_GLpYUEpetraDataPool.hpp"
00007 #include "rect2DMeshGenerator.hpp"
00008
00009 #include "Amesos_Klu.h"
00010 #include "EpetraExt_MatrixMatrix.h"
00011 #include "EpetraExt_Reindex_LinearProblem.h"
00012 #include "EpetraExt_Transpose_RowMatrix.h"
00013 #include "Epetra_BLAS.h"
00014 #include "Epetra_CrsMatrix.h"
00015 #include "Epetra_Export.h"
00016 #include "Epetra_FECrsMatrix.h"
00017 #include "Epetra_FEVector.h"
00018 #include "Epetra_Import.h"
00019 #include "Epetra_IntSerialDenseVector.h"
00020 #include "Epetra_LAPACK.h"
00021 #include "Epetra_Map.h"
00022 #include "Epetra_MultiVector.h"
00023 #include "Epetra_SerialDenseMatrix.h"
00024 #include "Epetra_Vector.h"
00025 #include "Teuchos_ParameterList.hpp"
00026 #include "Teuchos_RefCountPtr.hpp"
00027 #include "Teuchos_TestForException.hpp"
00028 #include "Teuchos_VerboseObject.hpp"
00029
00030 #ifdef HAVE_MPI
00031 # include "Epetra_MpiComm.h"
00032 # include "mpi.h"
00033 #else
00034 # include "Epetra_SerialComm.h"
00035 #endif
00036
00037
00038
00039 using namespace std;
00040
00041
00042
00043
00044
00045
00046
00047 namespace GLpApp {
00048
00049
00050
00051
00052
00053 const double GLp_pi = 3.14159265358979323846;
00054
00055 ostream& operator <<(ostream &, const Usr_Par &);
00056
00057 bool CrsMatrix2MATLAB(const Epetra_CrsMatrix &, ostream &);
00058
00059 bool Vector2MATLAB( const Epetra_Vector &, ostream &);
00060
00061 bool FEVector2MATLAB( const Epetra_FEVector &, ostream &);
00062
00063 int quadrature(const int, const int, Epetra_SerialDenseMatrix &,
00064 Epetra_SerialDenseVector &);
00065
00066 int meshreader(
00067 const Epetra_Comm &,
00068 Epetra_IntSerialDenseVector &,
00069 Epetra_SerialDenseMatrix &,
00070 Epetra_IntSerialDenseVector &,
00071 Epetra_SerialDenseMatrix &,
00072 Epetra_IntSerialDenseMatrix &,
00073 Epetra_IntSerialDenseMatrix &,
00074 const char geomFileBase[],
00075 const bool trace = true,
00076 const bool dumpAll = false
00077 );
00078
00079 int lassembly(const Epetra_SerialDenseMatrix &,
00080 const Epetra_SerialDenseVector &,
00081 const Epetra_SerialDenseMatrix &,
00082 const Epetra_SerialDenseVector &,
00083 const Epetra_SerialDenseVector &,
00084 const Usr_Par &,
00085 Epetra_SerialDenseMatrix &,
00086 Epetra_SerialDenseVector &);
00087
00088 int assemblyFECrs(const Epetra_Comm &,
00089 const Epetra_IntSerialDenseVector &,
00090 const Epetra_SerialDenseMatrix &,
00091 const Epetra_IntSerialDenseVector &,
00092 const Epetra_SerialDenseMatrix &,
00093 const Epetra_IntSerialDenseMatrix &,
00094 const Epetra_IntSerialDenseMatrix &,
00095 Teuchos::RefCountPtr<Epetra_FECrsMatrix> &,
00096 Teuchos::RefCountPtr<Epetra_FEVector> &);
00097
00098 int assemblyFECrs(const Epetra_Comm &,
00099 const Epetra_IntSerialDenseVector &,
00100 const Epetra_SerialDenseMatrix &,
00101 const Epetra_IntSerialDenseVector &,
00102 const Epetra_SerialDenseMatrix &,
00103 const Epetra_IntSerialDenseMatrix &,
00104 const Epetra_IntSerialDenseMatrix &,
00105 Teuchos::RefCountPtr<Epetra_FECrsMatrix> &,
00106 Teuchos::RefCountPtr<Epetra_FEVector> &,
00107 bool);
00108
00109 int assemble(const Epetra_Comm &,
00110 const Epetra_IntSerialDenseVector &,
00111 const Epetra_SerialDenseMatrix &,
00112 const Epetra_IntSerialDenseVector &,
00113 const Epetra_SerialDenseMatrix &,
00114 const Epetra_IntSerialDenseMatrix &,
00115 const Epetra_IntSerialDenseMatrix &,
00116 Teuchos::RefCountPtr<Epetra_FECrsMatrix> &,
00117 Teuchos::RefCountPtr<Epetra_FECrsMatrix> &,
00118 Teuchos::RefCountPtr<Epetra_FEVector> &);
00119
00120 int assemble_bdry(
00121 const Epetra_Comm &Comm
00122 ,const Epetra_IntSerialDenseVector &ipindx
00123 ,const Epetra_SerialDenseMatrix &ipcoords
00124 ,const Epetra_IntSerialDenseVector &pindx
00125 ,const Epetra_SerialDenseMatrix &pcoords
00126 ,const Epetra_IntSerialDenseMatrix &t
00127 ,const Epetra_IntSerialDenseMatrix &e
00128 ,Teuchos::RefCountPtr<Epetra_FECrsMatrix> *B
00129 ,Teuchos::RefCountPtr<Epetra_FECrsMatrix> *R
00130 );
00131
00132 int nonlinvec(const Epetra_Comm &,
00133 const Epetra_IntSerialDenseVector &,
00134 const Epetra_SerialDenseMatrix &,
00135 const Epetra_IntSerialDenseVector &,
00136 const Epetra_SerialDenseMatrix &,
00137 const Epetra_IntSerialDenseMatrix &,
00138 const Teuchos::RefCountPtr<const Epetra_MultiVector> &,
00139 Teuchos::RefCountPtr<Epetra_FEVector> &);
00140
00141 int nonlinjac(const Epetra_Comm &,
00142 const Epetra_IntSerialDenseVector &,
00143 const Epetra_SerialDenseMatrix &,
00144 const Epetra_IntSerialDenseVector &,
00145 const Epetra_SerialDenseMatrix &,
00146 const Epetra_IntSerialDenseMatrix &,
00147 const Teuchos::RefCountPtr<const Epetra_MultiVector> &,
00148 Teuchos::RefCountPtr<Epetra_FECrsMatrix> &);
00149
00150 int nonlinhessvec(const Epetra_Comm &,
00151 const Epetra_IntSerialDenseVector &,
00152 const Epetra_SerialDenseMatrix &,
00153 const Epetra_IntSerialDenseVector &,
00154 const Epetra_SerialDenseMatrix &,
00155 const Epetra_IntSerialDenseMatrix &,
00156 const Teuchos::RefCountPtr<const Epetra_MultiVector> &,
00157 const Teuchos::RefCountPtr<const Epetra_MultiVector> &,
00158 const Teuchos::RefCountPtr<const Epetra_MultiVector> &,
00159 Teuchos::RefCountPtr<Epetra_FEVector> &);
00160
00161 int compproduct(Epetra_SerialDenseVector &, double *, double *);
00162
00163 int compproduct(Epetra_SerialDenseVector &, double *, double *, double *);
00164
00165 double determinant(const Epetra_SerialDenseMatrix &);
00166
00167 int inverse(const Epetra_SerialDenseMatrix &, Epetra_SerialDenseMatrix &);
00168
00169 int quadrature(
00170 const int, const int, Epetra_SerialDenseMatrix &,
00171 Epetra_SerialDenseVector &);
00172
00173 void gpfctn(const Epetra_SerialDenseVector & v, Epetra_SerialDenseVector & gv);
00174
00175 void g2pfctn(const Epetra_SerialDenseVector & , Epetra_SerialDenseVector & );
00176
00177 void gfctn(const Epetra_SerialDenseVector & , Epetra_SerialDenseVector & );
00178
00179
00180
00181
00182
00183 GLpYUEpetraDataPool::GLpYUEpetraDataPool(
00184 Teuchos::RefCountPtr<const Epetra_Comm> const& commptr
00185 ,const double beta
00186 ,const double len_x
00187 ,const double len_y
00188 ,const int local_nx
00189 ,const int local_ny
00190 ,const char myfile[]
00191 ,const bool trace
00192 )
00193 :commptr_(commptr)
00194 ,beta_(beta)
00195 {
00196 ipcoords_ = Teuchos::rcp( new Epetra_SerialDenseMatrix() );
00197 ipindx_ = Teuchos::rcp( new Epetra_IntSerialDenseVector() );
00198 pcoords_ = Teuchos::rcp( new Epetra_SerialDenseMatrix() );
00199 pindx_ = Teuchos::rcp( new Epetra_IntSerialDenseVector() );
00200 t_ = Teuchos::rcp( new Epetra_IntSerialDenseMatrix() );
00201 e_ = Teuchos::rcp( new Epetra_IntSerialDenseMatrix() );
00202
00203 if( myfile && myfile[0] != '\0' ) {
00204 meshreader(*commptr_,*ipindx_,*ipcoords_,*pindx_,*pcoords_,*t_,*e_,myfile,trace);
00205 }
00206 else {
00207 rect2DMeshGenerator(
00208 commptr_->NumProc(),commptr_->MyPID()
00209 ,len_x,len_y,local_nx,local_ny,2
00210 ,&*ipindx_,&*ipcoords_,&*pindx_,&*pcoords_,&*t_,&*e_
00211 #ifdef GLPYUEPETRA_DATAPOOL_DUMP_ALL_MESH
00212 ,&*Teuchos::VerboseObjectBase::getDefaultOStream(),true
00213 #else
00214 ,NULL,false
00215 #endif
00216 );
00217 }
00218
00219
00220 assemble( *commptr, *ipindx_, *ipcoords_, *pindx_, *pcoords_, *t_, *e_, A_, H_, b_ );
00221 assemble_bdry( *commptr, *ipindx_, *ipcoords_, *pindx_, *pcoords_, *t_, *e_, &B_, &R_ );
00222
00223
00224 Epetra_Map standardmap(A_->DomainMap());
00225 q_ = Teuchos::rcp(new Epetra_FEVector(standardmap,1));
00226 int * qintvalues = new int[standardmap.NumMyElements()];
00227 double * qdvalues = new double[standardmap.NumMyElements()];
00228 standardmap.MyGlobalElements(qintvalues);
00229 for (int i = 0; i < standardmap.NumMyElements(); i++)
00230 qdvalues[i]=cos( GLp_pi* ((*ipcoords_)(i,0)) ) * cos( GLp_pi* ((*ipcoords_)(i,1)) );
00231 q_->ReplaceGlobalValues(standardmap.NumMyElements(), qintvalues, qdvalues);
00232 q_->GlobalAssemble();
00233 }
00234
00235 void GLpYUEpetraDataPool::computeAll( const GenSQP::Vector &x )
00236 {
00237
00238
00239 Teuchos::RefCountPtr<const Epetra_MultiVector> ey =
00240 (Teuchos::dyn_cast<GenSQP::YUEpetraVector>(const_cast<GenSQP::Vector&>(x))).getYVector();
00241
00242 computeNy(ey);
00243
00244 computeNpy(ey);
00245
00246 computeAugmat();
00247
00248 }
00249
00250 int GLpYUEpetraDataPool::solveAugsys(
00251 const Teuchos::RefCountPtr<const Epetra_MultiVector> & rhsy,
00252 const Teuchos::RefCountPtr<const Epetra_MultiVector> & rhsu,
00253 const Teuchos::RefCountPtr<const Epetra_MultiVector> & rhsp,
00254 const Teuchos::RefCountPtr<Epetra_MultiVector> & y,
00255 const Teuchos::RefCountPtr<Epetra_MultiVector> & u,
00256 const Teuchos::RefCountPtr<Epetra_MultiVector> & p,
00257 double tol )
00258 {
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272
00273
00274
00275
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325
00326
00327
00328
00329
00330
00331
00332
00333
00334
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344
00345
00346
00347
00348
00349
00350
00351
00352
00353
00354
00355
00356
00357
00358
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380
00381
00382
00383
00384
00385
00386
00387
00388
00389
00390
00391
00392
00393
00394
00395
00396
00397
00398
00399
00400
00401
00402
00403
00404
00405
00406
00407
00408
00409
00410
00411
00412
00413
00414
00415
00416
00417
00418
00419
00420
00421
00422
00423
00424
00425
00426
00427
00428
00429
00430
00431
00432
00433
00434
00435
00436
00437
00438
00439
00440
00441
00442
00443
00444
00445
00446
00447
00448
00449 TEST_FOR_EXCEPT(true);
00450 return 0;
00451 }
00452
00453 Teuchos::RefCountPtr<const Epetra_Comm> GLpYUEpetraDataPool::getCommPtr() { return commptr_; }
00454
00455 Teuchos::RefCountPtr<Epetra_FECrsMatrix> GLpYUEpetraDataPool::getA() { return A_; }
00456
00457 Teuchos::RefCountPtr<Epetra_FECrsMatrix> GLpYUEpetraDataPool::getB() { return B_; }
00458
00459 Teuchos::RefCountPtr<Epetra_FECrsMatrix> GLpYUEpetraDataPool::getH() { return H_; }
00460
00461 Teuchos::RefCountPtr<Epetra_FECrsMatrix> GLpYUEpetraDataPool::getR() { return R_; }
00462
00463 Teuchos::RefCountPtr<Epetra_CrsMatrix> GLpYUEpetraDataPool::getAugmat() { return Augmat_; }
00464
00465 Teuchos::RefCountPtr<Epetra_FECrsMatrix> GLpYUEpetraDataPool::getNpy() { return Npy_; }
00466
00467 Teuchos::RefCountPtr<Epetra_FEVector> GLpYUEpetraDataPool::getb() { return b_; }
00468
00469 Teuchos::RefCountPtr<Epetra_FEVector> GLpYUEpetraDataPool::getq() { return q_; }
00470
00471 Teuchos::RefCountPtr<Epetra_FEVector> GLpYUEpetraDataPool::getNy() { return Ny_; }
00472
00473 double GLpYUEpetraDataPool::getbeta() { return beta_; }
00474
00475 Teuchos::RefCountPtr<const Epetra_SerialDenseMatrix> GLpYUEpetraDataPool::getipcoords() { return ipcoords_; }
00476
00477 Teuchos::RefCountPtr<const Epetra_IntSerialDenseVector> GLpYUEpetraDataPool::getipindx() { return ipindx_; }
00478
00479 Teuchos::RefCountPtr<const Epetra_SerialDenseMatrix> GLpYUEpetraDataPool::getpcoords() { return pcoords_; }
00480
00481 Teuchos::RefCountPtr<const Epetra_IntSerialDenseVector> GLpYUEpetraDataPool::getpindx() { return pindx_; }
00482
00483 Teuchos::RefCountPtr<const Epetra_IntSerialDenseMatrix> GLpYUEpetraDataPool::gett() { return t_; }
00484
00485 Teuchos::RefCountPtr<const Epetra_IntSerialDenseMatrix> GLpYUEpetraDataPool::gete() { return e_; }
00486
00487
00488 void GLpYUEpetraDataPool::computeNy( const Teuchos::RefCountPtr<const Epetra_MultiVector> & y )
00489 {
00490 Epetra_Map overlapmap(-1, pindx_->M(), (int*)(pindx_)->A(), 1, *commptr_);
00491 Epetra_Map standardmap(A_->DomainMap());
00492 Teuchos::RefCountPtr<Epetra_MultiVector> yoverlap = Teuchos::rcp(new Epetra_MultiVector(overlapmap, 1));
00493 Epetra_Import Importer(overlapmap, standardmap);
00494 yoverlap->Import(*y, Importer, Insert);
00495 nonlinvec(*commptr_, *ipindx_, *ipcoords_, *pindx_, *pcoords_, *t_, yoverlap, Ny_);
00496 }
00497
00498
00499 void GLpYUEpetraDataPool::computeNpy( const Teuchos::RefCountPtr<const Epetra_MultiVector> & y )
00500 {
00501 Epetra_Map overlapmap(-1, pindx_->M(), (int*)(pindx_)->A(), 1, *commptr_);
00502 Epetra_Map standardmap(A_->DomainMap());
00503 Teuchos::RefCountPtr<Epetra_MultiVector> yoverlap = Teuchos::rcp(new Epetra_MultiVector(overlapmap, 1));
00504 Epetra_Import Importer(overlapmap, standardmap);
00505 yoverlap->Import(*y, Importer, Insert);
00506 nonlinjac(*commptr_, *ipindx_, *ipcoords_, *pindx_, *pcoords_, *t_, yoverlap, Npy_);
00507 }
00508
00509
00510 void GLpYUEpetraDataPool::computeAugmat()
00511 {
00512 Epetra_Map standardmap(A_->DomainMap());
00513 Epetra_Map bdryctrlmap(B_->DomainMap());
00514
00515 int indexBase = 1;
00516
00517 int numstates = standardmap.NumGlobalElements();
00518
00519 int nummystates = standardmap.NumMyElements();
00520 int nummycontrols = bdryctrlmap.NumMyElements();
00521
00522 Epetra_IntSerialDenseVector KKTmapindx(2*nummystates+nummycontrols);
00523
00524
00525 Epetra_IntSerialDenseVector states(nummystates);
00526 Epetra_IntSerialDenseVector controls(nummycontrols);
00527 standardmap.MyGlobalElements(states.Values());
00528 bdryctrlmap.MyGlobalElements(controls.Values());
00529 for (int i=0; i<nummystates; i++) {
00530 KKTmapindx(i) = states(i);
00531 KKTmapindx(nummystates+nummycontrols+i) = 2*numstates+states(i);
00532 }
00533 for (int i=0; i<nummycontrols; i++) {
00534 KKTmapindx(nummystates+i) = numstates+controls(i);
00535 }
00536 Epetra_Map KKTmap(-1, 2*nummystates+nummycontrols, KKTmapindx.Values(), indexBase, *(commptr_));
00537
00538
00539
00540
00541 Augmat_ = Teuchos::rcp(new Epetra_CrsMatrix(Copy, KKTmap, 0));
00542
00543 double one[1];
00544 one[0]=1.0;
00545 for (int i=0; i<nummystates+nummycontrols; i++) {
00546 Augmat_->InsertGlobalValues(KKTmapindx.Values()[i], 1, one, KKTmapindx.Values()+i);
00547 }
00548
00549 int checkentries=0;
00550 int nummyentries=0;
00551 Epetra_SerialDenseVector values(nummyentries);
00552 Epetra_IntSerialDenseVector indices(nummyentries);
00553
00554 for (int i=0; i<nummystates; i++) {
00555 nummyentries = A_->NumMyEntries(i);
00556 values.Resize(nummyentries);
00557 indices.Resize(nummyentries);
00558 A_->ExtractGlobalRowCopy(KKTmapindx.Values()[i], nummyentries, checkentries, values.Values(),
00559 indices.Values());
00560 if (nummyentries > 0)
00561 Augmat_->InsertGlobalValues(KKTmapindx.Values()[i]+2*numstates, nummyentries, values.Values(),
00562 indices.Values());
00563 nummyentries = Npy_->NumMyEntries(i);
00564 values.Resize(nummyentries);
00565 indices.Resize(nummyentries);
00566 Npy_->ExtractGlobalRowCopy(KKTmapindx.Values()[i], nummyentries, checkentries, values.Values(),
00567 indices.Values());
00568 if (nummyentries > 0)
00569 Augmat_->InsertGlobalValues(KKTmapindx.Values()[i]+2*numstates, nummyentries, values.Values(),
00570 indices.Values());
00571 }
00572
00573 for (int i=0; i<nummystates; i++) {
00574 nummyentries = B_->NumMyEntries(i);
00575 values.Resize(nummyentries);
00576 indices.Resize(nummyentries);
00577 B_->ExtractGlobalRowCopy(KKTmapindx.Values()[i], nummyentries, checkentries, values.Values(),
00578 indices.Values());
00579 for (int j=0; j<nummyentries; j++)
00580 indices[j] = indices[j]+numstates;
00581 if (nummyentries > 0)
00582 Augmat_->InsertGlobalValues(KKTmapindx.Values()[i]+2*numstates, nummyentries, values.Values(),
00583 indices.Values());
00584 }
00585
00586 bool MakeDataContiguous = false;
00587 EpetraExt::RowMatrix_Transpose transposer( MakeDataContiguous );
00588 Epetra_CrsMatrix & transA = dynamic_cast<Epetra_CrsMatrix&>(transposer(*A_));
00589 Epetra_CrsMatrix & transB = dynamic_cast<Epetra_CrsMatrix&>(transposer(*B_));
00590 Epetra_CrsMatrix & transNpy = dynamic_cast<Epetra_CrsMatrix&>(transposer(*Npy_));
00591
00592 for (int i=0; i<nummystates; i++) {
00593 nummyentries = transA.NumMyEntries(i);
00594 values.Resize(nummyentries);
00595 indices.Resize(nummyentries);
00596 transA.ExtractGlobalRowCopy(KKTmapindx.Values()[i], nummyentries, checkentries, values.Values(),
00597 indices.Values());
00598 for (int j=0; j<nummyentries; j++)
00599 indices[j] = indices[j]+2*numstates;
00600 if (nummyentries > 0)
00601 Augmat_->InsertGlobalValues(KKTmapindx.Values()[i], nummyentries, values.Values(),
00602 indices.Values());
00603 nummyentries = transNpy.NumMyEntries(i);
00604 values.Resize(nummyentries);
00605 indices.Resize(nummyentries);
00606 transNpy.ExtractGlobalRowCopy(KKTmapindx.Values()[i], nummyentries, checkentries, values.Values(),
00607 indices.Values());
00608 for (int j=0; j<nummyentries; j++)
00609 indices[j] = indices[j]+2*numstates;
00610 if (nummyentries > 0)
00611 Augmat_->InsertGlobalValues(KKTmapindx.Values()[i], nummyentries, values.Values(),
00612 indices.Values());
00613 }
00614
00615 for (int i=0; i<nummystates; i++) {
00616 nummyentries = transB.NumMyEntries(i);
00617 values.Resize(nummyentries);
00618 indices.Resize(nummyentries);
00619 transB.ExtractGlobalRowCopy(KKTmapindx.Values()[i], nummyentries, checkentries, values.Values(),
00620 indices.Values());
00621 for (int j=0; j<nummyentries; j++)
00622 indices[j] = indices[j]+2*numstates;
00623 int err = 0;
00624 if (nummyentries > 0)
00625 err = Augmat_->InsertGlobalValues(KKTmapindx.Values()[i]+numstates, nummyentries,
00626 values.Values(), indices.Values());
00627
00628 if (err < 0) {
00629 cout << "Insertion of entries failed:\n";
00630 cout << indices;
00631 cout << nummyentries << endl;
00632 cout << "at row: " << KKTmapindx.Values()[i]+numstates << endl << endl;
00633 }
00634 }
00635
00636 Augmat_->FillComplete(KKTmap, KKTmap);
00637
00638
00639 }
00640
00641 void GLpYUEpetraDataPool::PrintVec( const Teuchos::RefCountPtr<const Epetra_Vector> & x )
00642 {
00643 Vector2MATLAB(*x, cout);
00644 }
00645
00646 Usr_Par::Usr_Par() {
00647 Epetra_BLAS blas;
00648 Epetra_SerialDenseVector product(4);
00649
00650
00651 quadrature(2,3,Nodes,Weights);
00652
00653
00654
00655 N.Shape(Nodes.M(),3);
00656 for (int i=0; i<Nodes.M(); i++) {
00657 N(i,0) = 1.0 - Nodes(i,0) - Nodes(i,1);
00658 N(i,1) = Nodes(i,0);
00659 N(i,2) = Nodes(i,1);
00660 }
00661
00662 Nx1.Shape(Nodes.M(),3);
00663 for (int i=0; i<Nodes.M(); i++) {
00664 Nx1(i,0) = -1.0;
00665 Nx1(i,1) = 1.0;
00666 Nx1(i,2) = 0.0;
00667 }
00668
00669 Nx2.Shape(Nodes.M(),3);
00670 for (int i=0; i<Nodes.M(); i++) {
00671 Nx2(i,0) = -1.0;
00672 Nx2(i,1) = 0.0;
00673 Nx2(i,2) = 1.0;
00674 }
00675
00676
00677
00678 S1.Shape(3,3);
00679 S1(0,0)= 1.0; S1(0,1)=-1.0; S1(0,2)= 0.0;
00680 S1(1,0)=-1.0; S1(1,1)= 1.0; S1(1,2)= 0.0;
00681 S1(2,0)= 0.0; S1(2,1)= 0.0; S1(2,2)= 0.0;
00682 S2.Shape(3,3);
00683 S2(0,0)= 2.0; S2(0,1)=-1.0; S2(0,2)=-1.0;
00684 S2(1,0)=-1.0; S2(1,1)= 0.0; S2(1,2)= 1.0;
00685 S2(2,0)=-1.0; S2(2,1)= 1.0; S2(2,2)= 0.0;
00686 S3.Shape(3,3);
00687 S3(0,0)= 1.0; S3(0,1)= 0.0; S3(0,2)=-1.0;
00688 S3(1,0)= 0.0; S3(1,1)= 0.0; S3(1,2)= 0.0;
00689 S3(2,0)=-1.0; S3(2,1)= 0.0; S3(2,2)= 1.0;
00690
00691
00692 Nw.Size(3);
00693 Nw.Multiply('T', 'N', 1.0, N, Weights, 0.0);
00694
00695
00696 NNw.Shape(3,3);
00697 for (int i=0; i<3; i++) {
00698 for (int j=0; j<3; j++) {
00699 compproduct(product, N[i], N[j]);
00700 NNw(i,j) = blas.DOT(Weights.M(), Weights.A(), product.A());
00701 }
00702 }
00703
00704
00705 NNNw = new Epetra_SerialDenseMatrix[3];
00706 NNNw[0].Shape(3,3); NNNw[1].Shape(3,3); NNNw[2].Shape(3,3);
00707 for (int i=0; i<3; i++) {
00708 for (int j=0; j<3; j++) {
00709 for (int k=0; k<3; k++) {
00710 compproduct(product, N[i], N[j], N[k]);
00711 NNNw[k](i,j) = blas.DOT(Weights.M(), Weights.A(), product.A());
00712 }
00713 }
00714 }
00715
00716
00717 NdNdx1Nw = new Epetra_SerialDenseMatrix[3];
00718 NdNdx1Nw[0].Shape(3,3); NdNdx1Nw[1].Shape(3,3); NdNdx1Nw[2].Shape(3,3);
00719 for (int i=0; i<3; i++) {
00720 for (int j=0; j<3; j++) {
00721 for (int k=0; k<3; k++) {
00722 compproduct(product, N[i], Nx1[j], N[k]);
00723 NdNdx1Nw[k](i,j) = blas.DOT(Weights.M(), Weights.A(), product.A());
00724 }
00725 }
00726 }
00727
00728
00729 NdNdx2Nw = new Epetra_SerialDenseMatrix[3];
00730 NdNdx2Nw[0].Shape(3,3); NdNdx2Nw[1].Shape(3,3); NdNdx2Nw[2].Shape(3,3);
00731 for (int i=0; i<3; i++) {
00732 for (int j=0; j<3; j++) {
00733 for (int k=0; k<3; k++) {
00734 compproduct(product, N[i], Nx2[j], N[k]);
00735 NdNdx2Nw[k](i,j) = blas.DOT(Weights.M(), Weights.A(), product.A());
00736 }
00737 }
00738 }
00739
00740 }
00741
00742 void Usr_Par::Print(ostream& os) const {
00743 os << endl;
00744 os << "\n\nQuadrature nodes:";
00745 os << "\n-----------------";
00746 Nodes.Print(os);
00747 os << "\n\nQuadrature weights:";
00748 os << "\n-------------------\n";
00749 Weights.Print(os);
00750 os << "\n\nNodal basis functions:";
00751 os << "\n----------------------";
00752 N.Print(os);
00753 os << "\n\nIntegrals of combinations of partial derivatives:";
00754 os << "\n-------------------------------------------------";
00755 S1.Print(os);
00756 S2.Print(os);
00757 S3.Print(os);
00758 os << "\n\nIntegrals of basis functions:";
00759 os << "\n-----------------------------\n";
00760 Nw.Print(os);
00761 os << "\n\nIntegrals of products N(i)*N(j):";
00762 os << "\n--------------------------------\n";
00763 NNw.Print(os);
00764 os << "\n\nIntegrals of products N(i)*N(j)*N(k):";
00765 os << "\n-------------------------------------\n";
00766 NNNw[0].Print(os); NNNw[1].Print(os); NNNw[2].Print(os);
00767 os << "\n\nIntegrals of products N(i)*(dN(j)/dx1)*N(k):";
00768 os << "\n--------------------------------------------\n";
00769 NdNdx1Nw[0].Print(os); NdNdx1Nw[1].Print(os); NdNdx1Nw[2].Print(os);
00770 os << "\n\nIntegrals of products N(i)*(dN(j)/dx2)*N(k):";
00771 os << "\n--------------------------------------------\n";
00772 NdNdx2Nw[0].Print(os); NdNdx2Nw[1].Print(os); NdNdx2Nw[2].Print(os);
00773 }
00774
00775 ostream& operator <<(ostream & out, const Usr_Par & usr_par)
00776 {
00777 usr_par.Print(out);
00778 return out;
00779 }
00780
00781 }
00782
00783
00784
00785
00786
00787 int GLpApp::compproduct( Epetra_SerialDenseVector & product,
00788 double *first, double *second)
00789 {
00790 for (int i=0; i<product.M(); i++) {
00791 product[i] = first[i]*second[i];
00792 }
00793 return(0);
00794 }
00795
00796 int GLpApp::compproduct(Epetra_SerialDenseVector & product,
00797 double *first, double *second, double *third)
00798 {
00799 for (int i=0; i<product.M(); i++) {
00800 product[i] = first[i]*second[i]*third[i];
00801 }
00802 return(0);
00803 }
00804
00805
00806
00807
00808
00809
00810
00811
00812
00813
00814
00815
00816
00817
00818
00819
00820
00821
00822
00823
00824
00825
00826
00827
00828
00829
00830
00831
00832
00833
00834
00835
00836
00837
00838
00839
00840
00841
00842
00843
00844
00845
00846
00847
00848
00849 int GLpApp::assemble_bdry(
00850 const Epetra_Comm &Comm
00851 ,const Epetra_IntSerialDenseVector &ipindx
00852 ,const Epetra_SerialDenseMatrix &ipcoords
00853 ,const Epetra_IntSerialDenseVector &pindx
00854 ,const Epetra_SerialDenseMatrix &pcoords
00855 ,const Epetra_IntSerialDenseMatrix &t
00856 ,const Epetra_IntSerialDenseMatrix &e
00857 ,Teuchos::RefCountPtr<Epetra_FECrsMatrix> *B_out
00858 ,Teuchos::RefCountPtr<Epetra_FECrsMatrix> *R_out
00859 )
00860 {
00861
00862 using Teuchos::rcp;
00863
00864 #ifdef GLPAPP_SHOW_BOUNDARY_ASSEMBLY
00865 Teuchos::RefCountPtr<Teuchos::FancyOStream>
00866 out = Teuchos::VerboseObjectBase::getDefaultOStream();
00867 Teuchos::OSTab tab(out);
00868 *out << "\nEntering assemble_bdry(...) ...\n";
00869 #endif
00870
00871 int numLocElems = t.M();
00872 int numLocEdges = e.M();
00873
00874 int indexBase = 1;
00875
00876
00877
00878
00879 int * lastindx = 0;
00880 Epetra_IntSerialDenseVector BdryNode(2*numLocEdges);
00881 for (int j=0; j<numLocEdges; j++) {
00882 BdryNode[j] = e(j,0);
00883 BdryNode[j+numLocEdges] = e(j,1);
00884 }
00885 std::sort(BdryNode.Values(), BdryNode.Values()+2*numLocEdges);
00886 lastindx = std::unique(BdryNode.Values(), BdryNode.Values()+2*numLocEdges);
00887 const int numBdryNodes = lastindx - BdryNode.Values();
00888 BdryNode.Resize(numBdryNodes);
00889
00890
00891
00892
00893 Epetra_IntSerialDenseVector MyBdryNode(numBdryNodes);
00894 lastindx = std::set_intersection(
00895 BdryNode.Values(), BdryNode.Values()+numBdryNodes,
00896 ipindx.Values(), ipindx.Values()+ipindx.M(),
00897 MyBdryNode.Values()
00898 );
00899 const int numMyBdryNodes = lastindx - MyBdryNode.Values();
00900 MyBdryNode.Resize(numMyBdryNodes);
00901
00902
00903
00904
00905 Epetra_Map standardmap(-1, ipindx.M(), const_cast<int*>(ipindx.A()), indexBase, Comm);
00906 Epetra_Map overlapmap(-1, pindx.M(), const_cast<int*>(pindx.A()), indexBase, Comm);
00907 Epetra_Map mybdryctrlmap(-1, numMyBdryNodes, const_cast<int*>(MyBdryNode.A()), indexBase, Comm);
00908
00909
00910
00911
00912 #ifdef GLPAPP_SHOW_BOUNDARY_ASSEMBLY
00913 *out << "\nstandardmap:\n";
00914 standardmap.Print(Teuchos::OSTab(out).o());
00915 *out << "\nmybdryctrlmap:\n";
00916 mybdryctrlmap.Print(Teuchos::OSTab(out).o());
00917 #endif
00918
00919
00920
00921
00922 Teuchos::RefCountPtr<Epetra_FECrsMatrix>
00923 B = rcp(new Epetra_FECrsMatrix(Copy,standardmap,0)),
00924 R = rcp(new Epetra_FECrsMatrix(Copy,standardmap,0));
00925
00926
00927
00928 const int numNodesPerEdge = 2;
00929 Epetra_IntSerialDenseVector epetra_nodes(numNodesPerEdge);
00930
00931
00932
00933
00934
00935 Epetra_SerialDenseMatrix Bt(2,2);
00936 int err=0;
00937
00938 for( int i=0; i < numLocEdges; i++ ) {
00939
00940 const int
00941 global_id_0 = e(i,0),
00942 global_id_1 = e(i,1),
00943 local_id_0 = overlapmap.LID(global_id_0),
00944 local_id_1 = overlapmap.LID(global_id_1);
00945
00946 epetra_nodes(0) = global_id_0;
00947 epetra_nodes(1) = global_id_1;
00948
00949 const double
00950 x0 = pcoords(local_id_0,0),
00951 y0 = pcoords(local_id_0,1),
00952 x1 = pcoords(local_id_1,0),
00953 y1 = pcoords(local_id_1,1);
00954
00955 const double l = sqrt(pow(x0-x1,2) + pow(y0-y1,2));
00956
00957
00958 const double l_sixth = l * (1.0/6.0);
00959 Bt(0,0) = l_sixth * 2.0;
00960 Bt(0,1) = l_sixth * 1.0;
00961 Bt(1,0) = l_sixth * 1.0;
00962 Bt(1,1) = l_sixth * 2.0;
00963
00964 #ifdef GLPAPP_SHOW_BOUNDARY_ASSEMBLY
00965 *out
00966 << "\nEdge global nodes = ("<<global_id_0<<","<<global_id_1<<"),"
00967 << " local nodes = ("<<local_id_0<<","<<local_id_1<<"),"
00968 << " (x0,y0)(x1,y1)=("<<x0<<","<<y0<<")("<<x1<<","<<y1<<"),"
00969 << " Bt = ["<<Bt(0,0)<<","<<Bt(0,1)<<";"<<Bt(1,0)<<","<<Bt(1,1)<<"]\n";
00970 #endif
00971
00972 const int format = Epetra_FECrsMatrix::COLUMN_MAJOR;
00973 err = B->InsertGlobalValues(epetra_nodes,Bt,format);
00974 if (err<0) return(err);
00975 err = R->InsertGlobalValues(epetra_nodes,Bt,format);
00976 if (err<0) return(err);
00977
00978 }
00979
00980
00981
00982
00983
00984
00985
00986
00987
00988
00989
00990
00991
00992
00993
00994 err = B->GlobalAssemble(mybdryctrlmap,standardmap,true);
00995 if (err<0) return(err);
00996 err = R->GlobalAssemble(mybdryctrlmap,mybdryctrlmap,true);
00997 if (err<0) return(err);
00998
00999 if(B_out) *B_out = B;
01000 if(R_out) *R_out = R;
01001
01002 #ifdef GLPAPP_SHOW_BOUNDARY_ASSEMBLY
01003 *out << "\nB =\n";
01004 B->Print(Teuchos::OSTab(out).o());
01005 *out << "\nLeaving assemble_bdry(...) ...\n";
01006 #endif
01007
01008 return(0);
01009
01010 }
01011
01012
01013
01014
01015
01016
01017
01018
01019
01020
01021
01022
01023
01024
01025
01026
01027
01028
01029
01030
01031
01032
01033
01034
01035
01036
01037
01038
01039
01040
01041
01042
01043
01044
01045
01046
01047
01048
01049
01050
01051
01052
01053
01054
01055
01056
01057
01058
01059
01060
01061
01062
01063
01064
01065
01066
01067
01068 int GLpApp::assemble(const Epetra_Comm & Comm,
01069 const Epetra_IntSerialDenseVector & ipindx,
01070 const Epetra_SerialDenseMatrix & ipcoords,
01071 const Epetra_IntSerialDenseVector & pindx,
01072 const Epetra_SerialDenseMatrix & pcoords,
01073 const Epetra_IntSerialDenseMatrix & t,
01074 const Epetra_IntSerialDenseMatrix & e,
01075 RefCountPtr<Epetra_FECrsMatrix> & A,
01076 RefCountPtr<Epetra_FECrsMatrix> & M,
01077 RefCountPtr<Epetra_FEVector> & b)
01078 {
01079
01080 int myPID = Comm.MyPID();
01081 int numProcs = Comm.NumProc();
01082 Usr_Par usr_par;
01083
01084 int numLocElems = t.M();
01085 int numNodesPerElem = 3;
01086
01087 int indexBase = 1;
01088
01089 Epetra_Map standardmap(-1, ipindx.M(), (int*)ipindx.A(), indexBase, Comm);
01090 Epetra_Map overlapmap(-1, pindx.M(), (int*)pindx.A(), indexBase, Comm);
01091
01092 int* nodes = new int[numNodesPerElem];
01093 int i=0, j=0, err=0;
01094
01095 A = rcp(new Epetra_FECrsMatrix(Copy, standardmap, 0));
01096 M = rcp(new Epetra_FECrsMatrix(Copy, standardmap, 0));
01097 b = rcp(new Epetra_FEVector(standardmap,1));
01098
01099
01100 int format = Epetra_FECrsMatrix::COLUMN_MAJOR;
01101 Epetra_IntSerialDenseVector epetra_nodes(View, nodes, numNodesPerElem);
01102
01103
01104
01105 Epetra_SerialDenseMatrix At;
01106 Epetra_SerialDenseVector bt;
01107 Epetra_SerialDenseMatrix vertices(numNodesPerElem, pcoords.N());
01108
01109 Epetra_SerialDenseVector k(numNodesPerElem);
01110 for (i=0; i< numNodesPerElem; i++) k(i)=1.0;
01111 Epetra_SerialDenseMatrix c(numNodesPerElem,2);
01112 for (i=0; i< numNodesPerElem; i++) { c(i,0)=0.0; c(i,1)=0.0; }
01113 Epetra_SerialDenseVector r(numNodesPerElem);
01114 for (i=0; i< numNodesPerElem; i++) r(i)=0.0;
01115 Epetra_SerialDenseVector f(numNodesPerElem);
01116 for (i=0; i< numNodesPerElem; i++) f(i)=0.0;
01117 Epetra_SerialDenseVector g(2);
01118 g(0)=0.0; g(1)=0.0;
01119 Epetra_SerialDenseVector sigma(2);
01120 sigma(0)=0.0; sigma(1)=0.0;
01121 for(i=0; i<numLocElems; i++) {
01122 nodes[0] = t(i,0); nodes[1] = t(i,1); nodes[2] = t(i,2);
01123 for (j=0; j<numNodesPerElem; j++) {
01124
01125 vertices(j,0) = pcoords(overlapmap.LID(nodes[j]), 0);
01126 vertices(j,1) = pcoords(overlapmap.LID(nodes[j]), 1);
01127
01128 f(j) = cos(GLp_pi*vertices(j,0))*cos(GLp_pi*vertices(j,1)) *
01129 (2*GLp_pi*GLp_pi + pow(cos(GLp_pi*vertices(j,0)),2)*pow(cos(GLp_pi*vertices(j,1)),2) - 1.0);
01130 }
01131 lassembly(vertices, k, c, r, f, usr_par, At, bt);
01132 err = A->InsertGlobalValues(epetra_nodes, At, format);
01133 if (err<0) return(err);
01134 err = b->SumIntoGlobalValues(epetra_nodes, bt);
01135 if (err<0) return(err);
01136 }
01137
01138
01139
01140
01141 Epetra_IntSerialDenseMatrix neumann(e.M(),2);
01142 j = 0;
01143 for (i=0; i<e.M(); i++) {
01144 if (e(i,2)==2) {
01145 neumann(j,0) = e(i,0); neumann(j,1) = e(i,1);
01146 j++;
01147 }
01148 }
01149 neumann.Reshape(j,2);
01150
01151 if (neumann.M() != 0) {
01152 Epetra_BLAS blas;
01153 Epetra_SerialDenseMatrix quadnodes;
01154 Epetra_SerialDenseVector quadweights;
01155 Epetra_SerialDenseMatrix N;
01156 Epetra_SerialDenseMatrix NN;
01157 Epetra_SerialDenseVector product(2);
01158
01159
01160 quadrature(1, 2, quadnodes, quadweights);
01161
01162
01163
01164 N.Shape(quadnodes.M(),2);
01165 for (i=0; i<quadnodes.M(); i++) {
01166 N(i,0) = 1.0 - quadnodes(i,0);
01167 N(i,1) = quadnodes(i,0);
01168 }
01169
01170
01171 NN.Shape(2,2);
01172 for (i=0; i<2; i++) {
01173 for (j=0; j<2; j++) {
01174 compproduct(product, N[i], N[j]);
01175 NN(i,j) = blas.DOT(quadweights.M(), quadweights.A(), product.A());
01176 }
01177 }
01178
01179 Epetra_IntSerialDenseVector neumann_nodes(2);
01180 Epetra_SerialDenseVector neumann_add(2);
01181 double l;
01182 for (i=0; i<neumann.M(); i++) {
01183 neumann_nodes(0) = neumann(i,0); neumann_nodes(1) = neumann(i,1);
01184 neumann_add(0) = pcoords(overlapmap.LID(neumann_nodes(0)),0)
01185 -pcoords(overlapmap.LID(neumann_nodes(1)),0);
01186 neumann_add(1) = pcoords(overlapmap.LID(neumann_nodes(0)),1)
01187 -pcoords(overlapmap.LID(neumann_nodes(1)),1);
01188 l = blas.NRM2(neumann_add.M(), neumann_add.A());
01189 neumann_add.Multiply('N', 'N', 1.0, NN, g, 0.0);
01190 neumann_add.Scale(l);
01191 err = b->SumIntoGlobalValues(neumann_nodes, neumann_add);
01192 if (err<0) return(err);
01193 }
01194 }
01195
01196
01197 Epetra_IntSerialDenseMatrix robin(e.M(),2);
01198 j = 0;
01199 for (i=0; i<e.M(); i++) {
01200 if (e(i,2)==3) {
01201 robin(j,0) = e(i,0); robin(j,1) = e(i,1);
01202 j++;
01203 }
01204 }
01205 robin.Reshape(j,2);
01206
01207 if (robin.M() != 0) {
01208 Epetra_BLAS blas;
01209 Epetra_SerialDenseMatrix quadnodes;
01210 Epetra_SerialDenseVector quadweights;
01211 Epetra_SerialDenseMatrix N;
01212 Epetra_SerialDenseMatrix NN;
01213 Epetra_SerialDenseMatrix * NNN;
01214 Epetra_SerialDenseVector product(2);
01215
01216
01217 quadrature(1, 2, quadnodes, quadweights);
01218
01219
01220
01221 N.Shape(quadnodes.M(),2);
01222 for (i=0; i<quadnodes.M(); i++) {
01223 N(i,0) = 1.0 - quadnodes(i,0);
01224 N(i,1) = quadnodes(i,0);
01225 }
01226
01227
01228 NN.Shape(2,2);
01229 for (i=0; i<2; i++) {
01230 for (j=0; j<2; j++) {
01231 compproduct(product, N[i], N[j]);
01232 NN(i,j) = blas.DOT(quadweights.M(), quadweights.A(), product.A());
01233 }
01234 }
01235
01236
01237 NNN = new Epetra_SerialDenseMatrix[2];
01238 NNN[0].Shape(2,2); NNN[1].Shape(2,2);
01239 for (i=0; i<2; i++) {
01240 for (j=0; j<2; j++) {
01241 for (int k=0; k<2; k++) {
01242 compproduct(product, N[i], N[j], N[k]);
01243 NNN[k](i,j) = blas.DOT(quadweights.M(), quadweights.A(),
01244 product.A());
01245 }
01246 }
01247 }
01248
01249 Epetra_IntSerialDenseVector robin_nodes(2);
01250 Epetra_SerialDenseVector robin_b_add(2);
01251 Epetra_SerialDenseMatrix robin_A_add(2,2);
01252 double l;
01253 for (i=0; i<robin.M(); i++) {
01254 robin_nodes(0) = robin(i,0); robin_nodes(1) = robin(i,1);
01255
01256 robin_b_add(0) = pcoords(overlapmap.LID(robin_nodes(0)),0)
01257 -pcoords(overlapmap.LID(robin_nodes(1)),0);
01258 robin_b_add(1) = pcoords(overlapmap.LID(robin_nodes(0)),1)
01259 -pcoords(overlapmap.LID(robin_nodes(1)),1);
01260 l = blas.NRM2(robin_b_add.M(), robin_b_add.A());
01261 robin_b_add.Multiply('N', 'N', 1.0, NN, g, 0.0);
01262 robin_b_add.Scale(l);
01263 err = b->SumIntoGlobalValues(robin_nodes, robin_b_add);
01264 if (err<0) return(err);
01265
01266 NNN[0].Scale(sigma(0)); NNN[1].Scale(sigma(1));
01267 robin_A_add += NNN[0]; robin_A_add += NNN[1];
01268 robin_A_add.Scale(l);
01269 err = A->InsertGlobalValues(robin_nodes, robin_A_add, format);
01270 if (err<0) return(err);
01271 }
01272
01273 delete [] NNN;
01274 }
01275
01276
01277 Epetra_IntSerialDenseMatrix dirichlet(e.M(),2);
01278 j = 0;
01279 for (i=0; i<e.M(); i++) {
01280 if (e(i,2)==1) {
01281 dirichlet(j,0) = e(i,0); dirichlet(j,1) = e(i,1);
01282 j++;
01283 }
01284 }
01285 dirichlet.Reshape(j,2);
01286
01287
01288
01289
01290
01291
01292
01293
01294 Epetra_SerialDenseMatrix Mt;
01295
01296 for (i=0; i< numNodesPerElem; i++) k(i)=0.0;
01297 for (i=0; i< numNodesPerElem; i++) { c(i,0)=0.0; c(i,1)=0.0; }
01298 for (i=0; i< numNodesPerElem; i++) r(i)=1.0;
01299 for (i=0; i< numNodesPerElem; i++) f(i)=0.0;
01300 g(0)=0.0; g(1)=0.0;
01301 sigma(0)=0.0; sigma(1)=0.0;
01302 for(i=0; i<numLocElems; i++) {
01303 nodes[0] = t(i,0); nodes[1] = t(i,1); nodes[2] = t(i,2);
01304 for (j=0; j<numNodesPerElem; j++) {
01305 vertices(j,0) = pcoords(overlapmap.LID(nodes[j]), 0);
01306 vertices(j,1) = pcoords(overlapmap.LID(nodes[j]), 1);
01307 }
01308 lassembly(vertices, k, c, r, f, usr_par, Mt, bt);
01309 err = M->InsertGlobalValues(epetra_nodes, Mt, format);
01310 if (err<0) return(err);
01311 }
01312
01313
01314
01315
01316
01317
01318
01319 err = A->GlobalAssemble();
01320 if (err<0) return(err);
01321
01322 err = b->GlobalAssemble();
01323 if (err<0) return(err);
01324
01325 err = M->GlobalAssemble();
01326 if (err<0) return(err);
01327
01328 delete [] nodes;
01329
01330 return(0);
01331 }
01332
01333
01334
01335
01336
01337
01338
01339
01340
01341
01342
01343
01344
01345
01346
01347
01348
01349
01350
01351
01352
01353
01354
01355
01356
01357
01358
01359
01360
01361
01362
01363
01364
01365
01366
01367
01368
01369
01370
01371
01372
01373 int GLpApp::lassembly(const Epetra_SerialDenseMatrix & vertices,
01374 const Epetra_SerialDenseVector & k,
01375 const Epetra_SerialDenseMatrix & c,
01376 const Epetra_SerialDenseVector & r,
01377 const Epetra_SerialDenseVector & f,
01378 const Usr_Par & usr_par,
01379 Epetra_SerialDenseMatrix & At,
01380 Epetra_SerialDenseVector & bt)
01381 {
01382
01383 Epetra_SerialDenseMatrix B(2,2);
01384 Epetra_SerialDenseMatrix b(2,2);
01385 Epetra_SerialDenseMatrix BtB(2,2);
01386 Epetra_SerialDenseMatrix C(2,2);
01387 Epetra_SerialDenseMatrix M1(3,3);
01388 Epetra_SerialDenseMatrix M2(3,3);
01389 Epetra_SerialDenseMatrix M3(3,3);
01390 Epetra_SerialDenseMatrix tmp(3,3);
01391 double dB, adB;
01392 At.Shape(3,3);
01393 bt.Size(3);
01394
01395
01396 for(int i=0; i<2; i++) {
01397 B(i,0) = vertices(1,i)-vertices(0,i);
01398 B(i,1) = vertices(2,i)-vertices(0,i);
01399 }
01400 dB = determinant(B);
01401 adB = abs(dB);
01402
01403
01404 BtB.Multiply('T', 'N', 1.0, B, B, 0.0);
01405 inverse(BtB, C);
01406
01407 inverse(B, b); b.Scale(dB);
01408
01409
01410 tmp = usr_par.S1; tmp.Scale(C(0,0));
01411 M1 += tmp;
01412 tmp = usr_par.S2; tmp.Scale(C(0,1));
01413 M1 += tmp;
01414 tmp = usr_par.S3; tmp.Scale(C(1,1));
01415 M1 += tmp;
01416 M1.Scale( (k(0)*usr_par.Nw(0) + k(1)*usr_par.Nw(1) +
01417 k(2)*usr_par.Nw(2)) * adB );
01418
01419
01420 tmp = usr_par.NdNdx1Nw[0]; tmp.Scale(b(0,0)*c(0,0)+b(0,1)*c(0,1));
01421 M2 += tmp;
01422 tmp = usr_par.NdNdx1Nw[1]; tmp.Scale(b(0,0)*c(1,0)+b(0,1)*c(1,1));
01423 M2 += tmp;
01424 tmp = usr_par.NdNdx1Nw[2]; tmp.Scale(b(0,0)*c(2,0)+b(0,1)*c(2,1));
01425 M2 += tmp;
01426 tmp = usr_par.NdNdx2Nw[0]; tmp.Scale(b(1,0)*c(0,0)+b(1,1)*c(0,1));
01427 M2 += tmp;
01428 tmp = usr_par.NdNdx2Nw[1]; tmp.Scale(b(1,0)*c(1,0)+b(1,1)*c(1,1));
01429 M2 += tmp;
01430 tmp = usr_par.NdNdx2Nw[2]; tmp.Scale(b(1,0)*c(2,0)+b(1,1)*c(2,1));
01431 M2 += tmp;
01432 M2.Scale(adB/dB);
01433
01434
01435 tmp = usr_par.NNNw[0]; tmp.Scale(r(0));
01436 M3 += tmp;
01437 tmp = usr_par.NNNw[1]; tmp.Scale(r(1));
01438 M3 += tmp;
01439 tmp = usr_par.NNNw[2]; tmp.Scale(r(2));
01440 M3 += tmp;
01441 M3.Scale(adB);
01442
01443 At += M1;
01444 At += M2;
01445 At += M3;
01446
01447 bt.Multiply('T', 'N', adB, usr_par.NNw, f, 0.0);
01448
01449 return(0);
01450 }
01451
01452
01453
01454
01455
01456
01457
01458
01459 int GLpApp::inverse(const Epetra_SerialDenseMatrix & mat,
01460 Epetra_SerialDenseMatrix & inv)
01461 {
01462 Epetra_LAPACK lapack;
01463 int dim = mat.M();
01464 int info;
01465 Epetra_IntSerialDenseVector ipiv(dim);
01466 Epetra_SerialDenseVector work(2*dim);
01467
01468 inv.Shape(dim,dim);
01469 inv = mat;
01470
01471 lapack.GETRF(dim, dim, inv.A(), dim, ipiv.A(), &info);
01472 lapack.GETRI(dim, inv.A(), dim, ipiv.A(), work.A(), &dim, &info);
01473
01474 return(0);
01475 }
01476
01477
01478
01479
01480
01481
01482
01483 double GLpApp::determinant(const Epetra_SerialDenseMatrix & mat)
01484 {
01485
01486 Epetra_LAPACK lapack;
01487 double det;
01488 int swaps;
01489 int dim = mat.M();
01490 int info;
01491 Epetra_IntSerialDenseVector ipiv(dim);
01492
01493 Epetra_SerialDenseMatrix mymat(mat);
01494 lapack.GETRF(dim, dim, mymat.A(), dim, ipiv.A(), &info);
01495
01496 det = 1.0;
01497 for (int i=0; i<dim; i++) {
01498 det *= mymat(i,i);
01499 }
01500
01501 swaps = 0;
01502 for (int i=0; i<dim; i++) {
01503 if ((ipiv[i]-1) != i)
01504 swaps++;
01505 }
01506
01507 det *= pow((double)(-1.0),swaps);
01508
01509 return(det);
01510 }
01511
01512 int GLpApp::meshreader(const Epetra_Comm & Comm,
01513 Epetra_IntSerialDenseVector & ipindx,
01514 Epetra_SerialDenseMatrix & ipcoords,
01515 Epetra_IntSerialDenseVector & pindx,
01516 Epetra_SerialDenseMatrix & pcoords,
01517 Epetra_IntSerialDenseMatrix & t,
01518 Epetra_IntSerialDenseMatrix & e,
01519 const char geomFileBase[],
01520 const bool trace,
01521 const bool dumpAll
01522 )
01523 {
01524 int MyPID = Comm.MyPID();
01525
01526 const int FileNameSize = 120;
01527 char FileName[FileNameSize];
01528 TEST_FOR_EXCEPT(static_cast<int>(std::strlen(geomFileBase) + 5) > FileNameSize);
01529 sprintf(FileName, "%s.%03d", geomFileBase, MyPID);
01530
01531 {
01532 std::ifstream file_in(FileName);
01533 TEST_FOR_EXCEPTION(
01534 file_in.eof(), std::logic_error
01535 ,"Error, the file \""<<FileName<<"\" could not be opened for input!"
01536 );
01537 }
01538
01539 FILE* fid;
01540 fid = fopen(FileName, "r");
01541
01542 if(trace) printf("\nReading node info from %s ...\n", FileName);
01543 int numip = 0, numcp = 0;
01544 fscanf(fid, "%d %d", &numip, &numcp);
01545 const int nump = numip + numcp;
01546 if(trace) printf( "\nnumip = %d, numcp = %d, nump = %d\n", numip, numcp, nump );
01547 ipindx.Size(numip);
01548 ipcoords.Shape(numip, 2);
01549 pindx.Size(nump);
01550 pcoords.Shape(nump, 2);
01551 for (int i=0; i<numip; i++) {
01552 fscanf(fid,"%d %lf %lf %*d",&ipindx(i),&ipcoords(i,0),&ipcoords(i,1));
01553 if(trace&&dumpAll) printf("%d %lf %lf\n",ipindx(i),ipcoords(i,0),ipcoords(i,1));
01554 pindx(i) = ipindx(i);
01555 pcoords(i,0) = ipcoords(i,0); pcoords(i,1) = ipcoords(i,1);
01556 }
01557 for (int i=numip; i<nump; i++) {
01558 fscanf(fid,"%d %lf %lf %*d",&pindx(i),&pcoords(i,0),&pcoords(i,1));
01559 }
01560
01561 fscanf(fid,"%*[^\n]");
01562 fscanf(fid,"%*1[\n]");
01563
01564 fscanf(fid,"%*[^\n]");
01565 fscanf(fid,"%*1[\n]");
01566
01567 for (int i=0; i<nump; i++) {
01568 fscanf(fid,"%*[^\n]");
01569 fscanf(fid,"%*1[\n]");
01570 }
01571
01572 if(trace) printf("\nReading element info from %s ...\n", FileName);
01573 int numelems = 0;
01574 fscanf(fid, "%d", &numelems);
01575 if(trace) printf( "\nnumelems = %d\n", numelems );
01576 t.Shape(numelems,3);
01577 for (int i=0; i<numelems; i++) {
01578 fscanf(fid, "%d %d %d", &t(i,0), &t(i,1), &t(i,2));
01579 if(trace&&dumpAll) printf("%d %d %d\n", t(i,0), t(i,1), t(i,2));
01580 }
01581
01582 if(trace) printf("\nReading boundry edge info from %s ...\n", FileName);
01583 int numedges = 0;
01584 fscanf(fid,"%d",&numedges);
01585 if(trace) printf( "\nnumedges = %d\n", numedges );
01586 e.Shape(numedges,3);
01587 for (int i=0; i<numedges; i++) {
01588 fscanf(fid, "%d %d %d", &e(i,0), &e(i,1), &e(i,2));
01589 if(trace&&dumpAll) printf("%d %d %d\n", e(i,0), e(i,1), e(i,2));
01590 }
01591
01592 fclose(fid);
01593 if(trace) printf("Done reading mesh.\n");
01594
01595 return(0);
01596
01597 }
01598
01599
01600
01601
01602
01603
01604
01605
01606
01607
01608
01609
01610
01611
01612
01613
01614
01615
01616
01617
01618
01619
01620
01621
01622
01623
01624
01625
01626
01627
01628
01629
01630
01631
01632
01633
01634
01635 int GLpApp::nonlinhessvec(const Epetra_Comm & Comm,
01636 const Epetra_IntSerialDenseVector & ipindx,
01637 const Epetra_SerialDenseMatrix & ipcoords,
01638 const Epetra_IntSerialDenseVector & pindx,
01639 const Epetra_SerialDenseMatrix & pcoords,
01640 const Epetra_IntSerialDenseMatrix & t,
01641 const Teuchos::RefCountPtr<const Epetra_MultiVector> & y,
01642 const Teuchos::RefCountPtr<const Epetra_MultiVector> & s,
01643 const Teuchos::RefCountPtr<const Epetra_MultiVector> & lambda,
01644 Teuchos::RefCountPtr<Epetra_FEVector> & hessvec)
01645 {
01646
01647 int myPID = Comm.MyPID();
01648 int numProcs = Comm.NumProc();
01649
01650 int numLocNodes = pindx.M();
01651 int numMyLocNodes = ipindx.M();
01652 int numLocElems = t.M();
01653 int numNodesPerElem = 3;
01654
01655 int indexBase = 1;
01656
01657 Epetra_Map standardmap(-1, numMyLocNodes, (int*)ipindx.A(), indexBase, Comm);
01658 Epetra_Map overlapmap(-1, numLocNodes, (int*)pindx.A(), indexBase, Comm);
01659
01660 hessvec = rcp(new Epetra_FEVector(standardmap,1));
01661
01662 int* nodes = new int[numNodesPerElem];
01663 int i=0, j=0, err=0;
01664
01665
01666 Epetra_SerialDenseMatrix Nodes;
01667 Epetra_SerialDenseVector Weights;
01668 quadrature(2,3,Nodes,Weights);
01669 int numQuadPts = Nodes.M();
01670
01671
01672
01673 Epetra_SerialDenseMatrix N;
01674 N.Shape(numQuadPts,3);
01675 for (int i=0; i<numQuadPts; i++) {
01676 N(i,0) = 1.0 - Nodes(i,0) - Nodes(i,1);
01677 N(i,1) = Nodes(i,0);
01678 N(i,2) = Nodes(i,1);
01679 }
01680
01681
01682 Epetra_IntSerialDenseVector epetra_nodes(View, nodes, numNodesPerElem);
01683 Epetra_SerialDenseMatrix vertices(numNodesPerElem, pcoords.N());
01684
01685 Epetra_SerialDenseVector ly;
01686 Epetra_SerialDenseVector Nly;
01687 Epetra_SerialDenseVector ls;
01688 Epetra_SerialDenseVector Nls;
01689 Epetra_SerialDenseVector llambda;
01690 Epetra_SerialDenseVector Nllambda;
01691 Epetra_SerialDenseVector lgfctn;
01692 Epetra_SerialDenseVector lgfctnNi;
01693 Epetra_SerialDenseVector lgfctnNls;
01694 Epetra_SerialDenseVector lhessvec;
01695
01696 ly.Size(numNodesPerElem);
01697 Nly.Size(numQuadPts);
01698 ls.Size(numNodesPerElem);
01699 Nls.Size(numQuadPts);
01700 llambda.Size(numNodesPerElem);
01701 Nllambda.Size(numQuadPts);
01702 lgfctn.Size(numQuadPts);
01703 lgfctnNi.Size(numQuadPts);
01704 lgfctnNls.Size(numQuadPts);
01705 lhessvec.Size(numNodesPerElem);
01706
01707 Epetra_SerialDenseMatrix B(2,2);
01708 double adB;
01709
01710 for(i=0; i<numLocElems; i++) {
01711
01712 nodes[0] = t(i,0); nodes[1] = t(i,1); nodes[2] = t(i,2);
01713 for (j=0; j<numNodesPerElem; j++) {
01714 vertices(j,0) = pcoords(overlapmap.LID(nodes[j]), 0);
01715 vertices(j,1) = pcoords(overlapmap.LID(nodes[j]), 1);
01716 }
01717
01718
01719 for(int i=0; i<2; i++) {
01720 B(i,0) = vertices(1,i)-vertices(0,i);
01721 B(i,1) = vertices(2,i)-vertices(0,i);
01722 }
01723 adB = abs(determinant(B));
01724
01725
01726 for (j=0; j<numNodesPerElem; j++) {
01727 ly(j) = (*((*y)(0)))[overlapmap.LID(nodes[j])];
01728 ls(j) = (*((*s)(0)))[overlapmap.LID(nodes[j])];
01729 llambda(j) = (*((*lambda)(0)))[overlapmap.LID(nodes[j])];
01730 }
01731
01732 Nly.Multiply('N', 'N', 1.0, N, ly, 0.0);
01733 Nls.Multiply('N', 'N', 1.0, N, ls, 0.0);
01734 Nllambda.Multiply('N', 'N', 1.0, N, llambda, 0.0);
01735 g2pfctn(Nly, lgfctn);
01736
01737 for (int i=0; i<numNodesPerElem; i++) {
01738 compproduct(lgfctnNi, lgfctn.A(), N[i]);
01739 compproduct(lgfctnNls, lgfctnNi.A(), Nls.A(), Nllambda.A());
01740 lhessvec(i) = adB*lgfctnNls.Dot(Weights);
01741 }
01742
01743 err = hessvec->SumIntoGlobalValues(epetra_nodes, lhessvec);
01744 if (err<0) return(err);
01745 }
01746
01747
01748
01749 err = hessvec->GlobalAssemble();
01750 if (err<0) return(err);
01751
01752 delete [] nodes;
01753
01754 return(0);
01755 }
01756
01757
01758
01759
01760
01761 void GLpApp::g2pfctn(const Epetra_SerialDenseVector & v, Epetra_SerialDenseVector & gv) {
01762 for (int i=0; i<v.M(); i++) {
01763 gv(i) = 6.0*v(i);
01764 }
01765 }
01766
01767
01768
01769
01770
01771
01772
01773
01774
01775
01776
01777
01778
01779
01780
01781
01782
01783
01784
01785
01786
01787
01788
01789
01790
01791
01792
01793
01794
01795
01796
01797
01798
01799 int GLpApp::nonlinjac(const Epetra_Comm & Comm,
01800 const Epetra_IntSerialDenseVector & ipindx,
01801 const Epetra_SerialDenseMatrix & ipcoords,
01802 const Epetra_IntSerialDenseVector & pindx,
01803 const Epetra_SerialDenseMatrix & pcoords,
01804 const Epetra_IntSerialDenseMatrix & t,
01805 const Teuchos::RefCountPtr<const Epetra_MultiVector> & y,
01806 Teuchos::RefCountPtr<Epetra_FECrsMatrix> & Gp)
01807 {
01808
01809 int myPID = Comm.MyPID();
01810 int numProcs = Comm.NumProc();
01811
01812 int numLocNodes = pindx.M();
01813 int numMyLocNodes = ipindx.M();
01814 int numLocElems = t.M();
01815 int numNodesPerElem = 3;
01816
01817 int indexBase = 1;
01818
01819 Epetra_Map standardmap(-1, numMyLocNodes, (int*)ipindx.A(), indexBase, Comm);
01820 Epetra_Map overlapmap(-1, numLocNodes, (int*)pindx.A(), indexBase, Comm);
01821
01822 int format = Epetra_FECrsMatrix::COLUMN_MAJOR;
01823 Gp = rcp(new Epetra_FECrsMatrix(Copy, standardmap, 0));
01824
01825 int* nodes = new int[numNodesPerElem];
01826 int i=0, j=0, err=0;
01827
01828
01829 Epetra_SerialDenseMatrix Nodes;
01830 Epetra_SerialDenseVector Weights;
01831 quadrature(2,3,Nodes,Weights);
01832 int numQuadPts = Nodes.M();
01833
01834
01835
01836 Epetra_SerialDenseMatrix N;
01837 N.Shape(numQuadPts,3);
01838 for (int i=0; i<numQuadPts; i++) {
01839 N(i,0) = 1.0 - Nodes(i,0) - Nodes(i,1);
01840 N(i,1) = Nodes(i,0);
01841 N(i,2) = Nodes(i,1);
01842 }
01843
01844
01845 Epetra_IntSerialDenseVector epetra_nodes(View, nodes, numNodesPerElem);
01846 Epetra_SerialDenseMatrix vertices(numNodesPerElem, pcoords.N());
01847
01848 Epetra_SerialDenseVector ly;
01849 Epetra_SerialDenseVector Nly;
01850 Epetra_SerialDenseVector lgfctn;
01851 Epetra_SerialDenseVector lgfctnNiNj;
01852 Epetra_SerialDenseMatrix lGp;
01853
01854 ly.Size(numNodesPerElem);
01855 Nly.Size(numQuadPts);
01856 lgfctn.Size(numQuadPts);
01857 lgfctnNiNj.Size(numQuadPts);
01858 lGp.Shape(numNodesPerElem, numNodesPerElem);
01859
01860 Epetra_SerialDenseMatrix B(2,2);
01861 double adB;
01862
01863 for(i=0; i<numLocElems; i++) {
01864
01865 nodes[0] = t(i,0); nodes[1] = t(i,1); nodes[2] = t(i,2);
01866 for (j=0; j<numNodesPerElem; j++) {
01867 vertices(j,0) = pcoords(overlapmap.LID(nodes[j]), 0);
01868 vertices(j,1) = pcoords(overlapmap.LID(nodes[j]), 1);
01869 }
01870
01871
01872 for(int i=0; i<2; i++) {
01873 B(i,0) = vertices(1,i)-vertices(0,i);
01874 B(i,1) = vertices(2,i)-vertices(0,i);
01875 }
01876 adB = abs(determinant(B));
01877
01878
01879 for (j=0; j<numNodesPerElem; j++) {
01880 ly(j) = (*((*y)(0)))[overlapmap.LID(nodes[j])];
01881 }
01882
01883 Nly.Multiply('N', 'N', 1.0, N, ly, 0.0);
01884 gpfctn(Nly, lgfctn);
01885
01886 for (int i=0; i<numNodesPerElem; i++) {
01887 for (int j=0; j<numNodesPerElem; j++) {
01888 compproduct(lgfctnNiNj, lgfctn.A(), N[i], N[j]);
01889 lGp(i,j) = adB*lgfctnNiNj.Dot(Weights);
01890 }
01891 }
01892
01893 err = Gp->InsertGlobalValues(epetra_nodes, lGp, format);
01894 if (err<0) return(err);
01895 }
01896
01897
01898
01899 err = Gp->GlobalAssemble();
01900 if (err<0) return(err);
01901
01902 delete [] nodes;
01903
01904 return(0);
01905 }
01906
01907
01908
01909
01910
01911 void GLpApp::gpfctn(const Epetra_SerialDenseVector & v, Epetra_SerialDenseVector & gv) {
01912 for (int i=0; i<v.M(); i++) {
01913 gv(i) = 3.0*pow(v(i),2)-1.0;
01914 }
01915 }
01916
01917
01918
01919
01920
01921
01922
01923
01924
01925
01926
01927
01928
01929
01930
01931
01932
01933
01934
01935
01936
01937
01938
01939
01940
01941
01942
01943
01944
01945
01946
01947
01948
01949 int GLpApp::nonlinvec(const Epetra_Comm & Comm,
01950 const Epetra_IntSerialDenseVector & ipindx,
01951 const Epetra_SerialDenseMatrix & ipcoords,
01952 const Epetra_IntSerialDenseVector & pindx,
01953 const Epetra_SerialDenseMatrix & pcoords,
01954 const Epetra_IntSerialDenseMatrix & t,
01955 const Teuchos::RefCountPtr<const Epetra_MultiVector> & y,
01956 Teuchos::RefCountPtr<Epetra_FEVector> & g)
01957 {
01958
01959 int myPID = Comm.MyPID();
01960 int numProcs = Comm.NumProc();
01961
01962 int numLocNodes = pindx.M();
01963 int numMyLocNodes = ipindx.M();
01964 int numLocElems = t.M();
01965 int numNodesPerElem = 3;
01966
01967 int indexBase = 1;
01968
01969 Epetra_Map standardmap(-1, numMyLocNodes, (int*)ipindx.A(), indexBase, Comm);
01970 Epetra_Map overlapmap(-1, numLocNodes, (int*)pindx.A(), indexBase, Comm);
01971
01972 g = rcp(new Epetra_FEVector(standardmap,1));
01973
01974 int* nodes = new int[numNodesPerElem];
01975 int i=0, j=0, err=0;
01976
01977
01978 Epetra_SerialDenseMatrix Nodes;
01979 Epetra_SerialDenseVector Weights;
01980 quadrature(2,3,Nodes,Weights);
01981 int numQuadPts = Nodes.M();
01982
01983
01984
01985 Epetra_SerialDenseMatrix N;
01986 N.Shape(numQuadPts,3);
01987 for (int i=0; i<numQuadPts; i++) {
01988 N(i,0) = 1.0 - Nodes(i,0) - Nodes(i,1);
01989 N(i,1) = Nodes(i,0);
01990 N(i,2) = Nodes(i,1);
01991 }
01992
01993
01994 Epetra_IntSerialDenseVector epetra_nodes(View, nodes, numNodesPerElem);
01995 Epetra_SerialDenseMatrix vertices(numNodesPerElem, pcoords.N());
01996
01997 Epetra_SerialDenseVector ly;
01998 Epetra_SerialDenseVector Nly;
01999 Epetra_SerialDenseVector lgfctn;
02000 Epetra_SerialDenseVector lgfctnNi;
02001 Epetra_SerialDenseVector lg;
02002
02003 ly.Size(numNodesPerElem);
02004 Nly.Size(numQuadPts);
02005 lgfctn.Size(numQuadPts);
02006 lgfctnNi.Size(numQuadPts);
02007 lg.Size(numNodesPerElem);
02008
02009 Epetra_SerialDenseMatrix B(2,2);
02010 double adB;
02011
02012 for(i=0; i<numLocElems; i++) {
02013
02014 nodes[0] = t(i,0); nodes[1] = t(i,1); nodes[2] = t(i,2);
02015 for (j=0; j<numNodesPerElem; j++) {
02016 vertices(j,0) = pcoords(overlapmap.LID(nodes[j]), 0);
02017 vertices(j,1) = pcoords(overlapmap.LID(nodes[j]), 1);
02018 }
02019
02020
02021 for(int i=0; i<2; i++) {
02022 B(i,0) = vertices(1,i)-vertices(0,i);
02023 B(i,1) = vertices(2,i)-vertices(0,i);
02024 }
02025 adB = abs(determinant(B));
02026
02027
02028 for (j=0; j<numNodesPerElem; j++) {
02029 ly(j) = (*((*y)(0)))[overlapmap.LID(nodes[j])];
02030 }
02031
02032 Nly.Multiply('N', 'N', 1.0, N, ly, 0.0);
02033 gfctn(Nly, lgfctn);
02034
02035 for (int i=0; i<numNodesPerElem; i++) {
02036 compproduct(lgfctnNi, lgfctn.A(), N[i]);
02037 lg(i) = adB*lgfctnNi.Dot(Weights);
02038 }
02039
02040 err = g->SumIntoGlobalValues(epetra_nodes, lg);
02041 if (err<0) return(err);
02042 }
02043
02044
02045
02046 err = g->GlobalAssemble();
02047 if (err<0) return(err);
02048
02049 delete [] nodes;
02050
02051 return(0);
02052 }
02053
02054
02055
02056
02057
02058
02059 void GLpApp::gfctn(const Epetra_SerialDenseVector & v, Epetra_SerialDenseVector & gv) {
02060 for (int i=0; i<v.M(); i++) {
02061 gv(i) = pow(v(i),3)-v(i);
02062 }
02063 }
02064
02065
02066
02067
02068
02069
02070
02071
02072
02073
02074
02075
02076
02077
02078
02079
02080
02081
02082
02083
02084
02085
02086 bool GLpApp::CrsMatrix2MATLAB(const Epetra_CrsMatrix & A, ostream & outfile)
02087 {
02088
02089 int MyPID = A.Comm().MyPID();
02090 int NumProc = A.Comm().NumProc();
02091
02092
02093 if( A.IndicesAreLocal() == false ) {
02094 if( MyPID == 0 ) {
02095 cerr << "ERROR in "<< __FILE__ << ", line " << __LINE__ << endl;
02096 cerr << "Function CrsMatrix2MATLAB accepts\n";
02097 cerr << "transformed matrices ONLY. Please call A.TransformToLoca()\n";
02098 cerr << "on your matrix A to that purpose.\n";
02099 cerr << "Now returning...\n";
02100 }
02101 return false;
02102 }
02103
02104 int NumMyRows = A.NumMyRows();
02105 int NumNzRow;
02106 int NumEntries;
02107 int NumGlobalRows;
02108 int GlobalRow;
02109 int NumGlobalNonzeros;
02110
02111 NumGlobalRows = A.NumGlobalRows();
02112 NumGlobalNonzeros = A.NumGlobalNonzeros();
02113
02114
02115
02116 int IndexBase = A.IndexBase();
02117 if( IndexBase == 0 )
02118 IndexBase = 1;
02119 else if ( IndexBase == 1)
02120 IndexBase = 0;
02121
02122
02123
02124 if( MyPID==0 ) {
02125 outfile << "A = spalloc(";
02126 outfile << NumGlobalRows << ',' << NumGlobalRows;
02127 outfile << ',' << NumGlobalNonzeros << ");\n";
02128 }
02129
02130 for( int Proc=0 ; Proc<NumProc ; ++Proc ) {
02131 A.Comm().Barrier();
02132 if( MyPID == Proc ) {
02133
02134 outfile << "\n\n% On proc " << Proc << ": ";
02135 outfile << NumMyRows << " rows and ";
02136 outfile << A.NumMyNonzeros() << " nonzeros\n";
02137
02138
02139 for( int MyRow=0 ; MyRow<NumMyRows ; ++MyRow ) {
02140
02141 GlobalRow = A.GRID(MyRow);
02142
02143 NumNzRow = A.NumMyEntries(MyRow);
02144 double *Values = new double[NumNzRow];
02145 int *Indices = new int[NumNzRow];
02146
02147 A.ExtractMyRowCopy(MyRow, NumNzRow,
02148 NumEntries, Values, Indices);
02149
02150 for( int j=0 ; j<NumEntries ; ++j ) {
02151 outfile << "A(" << GlobalRow + IndexBase
02152 << "," << A.GCID(Indices[j]) + IndexBase
02153 << ") = " << Values[j] << ";\n";
02154 }
02155
02156 delete Values;
02157 delete Indices;
02158 }
02159
02160 }
02161 A.Comm().Barrier();
02162 }
02163
02164 return true;
02165
02166 }
02167
02168
02169
02170
02171
02172
02173
02174
02175
02176
02177
02178
02179
02180
02181
02182
02183
02184
02185
02186
02187
02188 bool GLpApp::Vector2MATLAB( const Epetra_Vector & v, ostream & outfile)
02189 {
02190
02191 int MyPID = v.Comm().MyPID();
02192 int NumProc = v.Comm().NumProc();
02193 int MyLength = v.MyLength();
02194 int GlobalLength = v.GlobalLength();
02195
02196
02197
02198
02199
02200 if( MyPID == 0 ) outfile << "v = zeros(" << GlobalLength << ",1)\n";
02201
02202 int NumMyElements = v.Map().NumMyElements();
02203
02204 int * MyGlobalElements = v.Map().MyGlobalElements( );
02205
02206 int Row;
02207
02208 int IndexBase = v.Map().IndexBase();
02209 if( IndexBase == 0 )
02210 IndexBase = 1;
02211 else if ( IndexBase == 1)
02212 IndexBase = 0;
02213
02214 for( int Proc=0 ; Proc<NumProc ; ++Proc ) {
02215 v.Comm().Barrier();
02216 if( MyPID == Proc ) {
02217
02218 outfile << "% On proc " << Proc << ": ";
02219 outfile << MyLength << " rows of ";
02220 outfile << GlobalLength << " elements\n";
02221
02222 for( Row=0 ; Row<MyLength ; ++Row ) {
02223 outfile << "v(" << MyGlobalElements[Row] + IndexBase
02224 << ") = " << v[Row] << ";\n";
02225 }
02226
02227 }
02228
02229 v.Comm().Barrier();
02230 }
02231
02232 return true;
02233
02234 }
02235
02236
02237
02238
02239
02240
02241
02242
02243
02244
02245
02246
02247
02248
02249
02250
02251
02252
02253
02254
02255
02256 bool GLpApp::FEVector2MATLAB( const Epetra_FEVector & v, ostream & outfile)
02257 {
02258
02259 int MyPID = v.Comm().MyPID();
02260 int NumProc = v.Comm().NumProc();
02261 int MyLength = v.MyLength();
02262 int GlobalLength = v.GlobalLength();
02263
02264
02265
02266
02267
02268 if( MyPID == 0 ) outfile << "v = zeros(" << GlobalLength << ",1)\n";
02269
02270 int NumMyElements = v.Map().NumMyElements();
02271
02272 int * MyGlobalElements = v.Map().MyGlobalElements( );
02273
02274 int Row;
02275
02276 int IndexBase = v.Map().IndexBase();
02277 if( IndexBase == 0 )
02278 IndexBase = 1;
02279 else if ( IndexBase == 1)
02280 IndexBase = 0;
02281
02282 for( int Proc=0 ; Proc<NumProc ; ++Proc ) {
02283 v.Comm().Barrier();
02284 if( MyPID == Proc ) {
02285
02286 outfile << "% On proc " << Proc << ": ";
02287 outfile << MyLength << " rows of ";
02288 outfile << GlobalLength << " elements\n";
02289
02290 for( Row=0 ; Row<MyLength ; ++Row ) {
02291 outfile << "v(" << MyGlobalElements[Row] + IndexBase
02292 << ") = " << v[0][Row] << ";\n";
02293 }
02294
02295 }
02296
02297 v.Comm().Barrier();
02298 }
02299
02300 return true;
02301
02302 }
02303
02304
02305
02306
02307
02308
02309
02310
02311
02312
02313
02314
02315
02316
02317 int GLpApp::quadrature(const int dim, const int order,
02318 Epetra_SerialDenseMatrix & nodes,
02319 Epetra_SerialDenseVector & weights)
02320 {
02321
02322 if (dim == 1) {
02323
02324
02325
02326 if (order == 1) {
02327 nodes.Shape(1,1);
02328 nodes(0,0) = 0.5;
02329 weights.Size(1);
02330 weights(0) = 1.0;
02331 }
02332 else if (order == 2) {
02333 nodes.Shape(2,1);
02334 nodes(0,0) = (1.0-1.0/sqrt(3.0))/2.0;
02335 nodes(1,0) = (1.0+1.0/sqrt(3.0))/2.0;
02336 weights.Size(2);
02337 weights(0) = 0.5;
02338 weights(1) = 0.5;
02339 }
02340 else if (order == 3) {
02341 nodes.Shape(3,1);
02342 nodes(0,0) = (1.0-sqrt(3.0/5.0))/2.0;
02343 nodes(1,0) = 0.5;
02344 nodes(2,0) = (1.0+sqrt(3.0/5.0))/2.0;
02345 weights.Size(3);
02346 weights(0) = 5.0/18.0;
02347 weights(1) = 4.0/9.0;
02348 weights(2) = 5.0/18.0;
02349 }
02350 else {
02351 cout << "Quadrature for dim = " << dim << " and order = ";
02352 cout << order << " not available.\n";
02353 exit(-1);
02354 }
02355
02356 }
02357 else if (dim == 2) {
02358
02359
02360
02361
02362 if (order == 1) {
02363 nodes.Shape(1,2);
02364 nodes(0,0) = 1.0/3.0; nodes (0,1) = 1.0/3.0;
02365 weights.Size(1);
02366 weights(0) = 0.5;
02367 }
02368 else if (order == 2) {
02369 nodes.Shape(3,2);
02370 nodes(0,0) = 1.0/6.0; nodes (0,1) = 1.0/6.0;
02371 nodes(1,0) = 2.0/3.0; nodes (1,1) = 1.0/6.0;
02372 nodes(2,0) = 1.0/6.0; nodes (2,1) = 2.0/3.0;
02373 weights.Size(3);
02374 weights(0) = 1.0/6.0;
02375 weights(1) = 1.0/6.0;
02376 weights(2) = 1.0/6.0;
02377 }
02378 else if (order == 3) {
02379 nodes.Shape(4,2);
02380 nodes(0,0) = 1.0/3.0; nodes (0,1) = 1.0/3.0;
02381 nodes(1,0) = 3.0/5.0; nodes (1,1) = 1.0/5.0;
02382 nodes(2,0) = 1.0/5.0; nodes (2,1) = 3.0/5.0;
02383 nodes(3,0) = 1.0/5.0; nodes (3,1) = 1.0/5.0;
02384 weights.Size(4);
02385 weights(0) = -9.0/32.0;
02386 weights(1) = 25.0/96.0;
02387 weights(2) = 25.0/96.0;
02388 weights(3) = 25.0/96.0;
02389 }
02390 else {
02391 cout << "Quadrature for dim = " << dim << " and order = ";
02392 cout << order << " not available.\n";
02393 exit(-1);
02394 }
02395
02396 }
02397 else {
02398 cout << "Quadrature for dim = " << dim << " not available.\n";
02399 exit(-1);
02400 }
02401
02402 return(0);
02403 }