/*----------------------------------------------------------------------* | project master nodes onto slave element (private) mwgee 10/05| *----------------------------------------------------------------------*/ bool MOERTEL::Overlap::build_mxi() { // project the master segment's nodes onto the slave segment MOERTEL::Node** mnode = mseg_.Nodes(); MOERTEL::Projector projector(inter_.IsOneDimensional(),OutLevel()); double gap; for (int i=0; i<mseg_.Nnode(); ++i) { // project node i onto sseg projector.ProjectNodetoSegment_SegmentNormal(*mnode[i],sseg_,mxi_[i],gap); /* mxi_[i] output. mxi_[i][0] is the \xi coordinate of the projection in sseg_, and mxi_[i][1] is the \eta coordinate of the projection in sseg_. */ #if 0 // check whether i is inside sseg if (mxi_[i][0]<=1. && mxi_[i][1]<=abs(1.-mxi_[i][0]) && mxi_[i][0]>=0. && mxi_[i][1]>=0.) { std::cout << "OVERLAP: master node in: " << mnode[i]->Id() << endl; } #endif } havemxi_ = true; return true; }
/*----------------------------------------------------------------------* | add point (private) mwgee 10/05| *----------------------------------------------------------------------*/ bool MOERTEL::Overlap::AddPointtoPolygon(std::map<int,Teuchos::RCP<MOERTEL::Point> >& p, const int id,const double* P) { Teuchos::RCP<MOERTEL::Point> point = Teuchos::rcp(new MOERTEL::Point(id,P,OutLevel())); p.insert(std::pair<int,Teuchos::RCP<MOERTEL::Point> >(id,point)); return true; }
/*----------------------------------------------------------------------* | mwgee 08/05| *----------------------------------------------------------------------*/ bool MOERTEL::Projector::ProjectNodetoSegment_SegmentOrthogonal(MOERTEL::Node& node, MOERTEL::Segment& seg, double xi[], double &gap) { #if 0 std::cout << "----- Projector: Node " << node.Id() << " Segment " << seg.Id() << std::endl; #endif if (IsTwoDimensional()) { // we do a newton iteration for the projection coordinates xi // set starting value to the middle of the segment double eta = 0.0; int i = 0; double F=0.0,dF=0.0,deta=0.0; for (i=0; i<10; ++i) { F = evaluate_F_2D_SegmentOrthogonal(node,seg,eta,gap); if (abs(F) < 1.0e-10) break; dF = evaluate_gradF_2D_SegmentOrthogonal(node,seg,eta); deta = (-F)/dF; eta += deta; } if (abs(F)>1.0e-10) { if (OutLevel()>3) std::cout << "MOERTEL: ***WRN*** MOERTEL::Projector::ProjectNodetoSegment_SegmentOrthogonal:\n" << "MOERTEL: ***WRN*** Newton iteration failed to converge\n" << "MOERTEL: ***WRN*** #iterations = " << i << std::endl << "MOERTEL: ***WRN*** F(eta) = " << F << " gradF(eta) = " << dF << " eta = " << eta << " delta(eta) = " << deta << "\n" << "MOERTEL: ***WRN*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; } #if 0 std::cout << "#iterations = " << i << " F = " << F << " eta = " << eta << std::endl; #endif xi[0] = eta; return true; } else { std::stringstream oss; oss << "***ERR*** MOERTEL::Projector::ProjectNodetoSegment_SegmentOrthogonal:\n" << "***ERR*** 3D projection not yet impl.\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw ReportError(oss); } return true; }
/*----------------------------------------------------------------------* | project slave nodes onto master element (private) mwgee 11/05| *----------------------------------------------------------------------*/ bool MOERTEL::Overlap::build_sxim() { // project the slave segment's nodes onto the master segment int nsnode = sseg_.Nnode(); MOERTEL::Node** snode = sseg_.Nodes(); MOERTEL::Projector projector(inter_.IsOneDimensional(),OutLevel()); double gap; for (int i=0; i<nsnode; ++i) { // project node i onto sseg projector.ProjectNodetoSegment_NodalNormal(*snode[i],mseg_,sxim_[i],gap); #if 0 // check whether i is inside sseg if (sxim_[i][0]<=1. && sxim_[i][1]<=abs(1.-sxim_[i][0]) && sxim_[i][0]>=0. && sxim_[i][1]>=0.) { std::cout << "OVERLAP: slave node in: " << snode[i]->Id() << endl; } #endif } havesxim_ = true; return true; }
/*----------------------------------------------------------------------* | add point (private) mwgee 10/05| *----------------------------------------------------------------------*/ bool MOERTEL::Overlap::AddPointtoPolygon(const int id,const double* P) { // check whether this point is already in there std::map<int,Teuchos::RCP<MOERTEL::Point> >::iterator curr = p_.find(id); // it's there if (curr != p_.end()) { curr->second->SetXi(P); } else { // std::cout << "OVERLAP Clip_AddPointtoPolygon: added point " << id // << " xi=" << P[0] << "/" << P[1] << std::endl; Teuchos::RCP<MOERTEL::Point> p = Teuchos::rcp(new MOERTEL::Point(id,P,OutLevel())); p_.insert(std::pair<int,Teuchos::RCP<MOERTEL::Point> >(id,p)); } return true; }
/*----------------------------------------------------------------------* | finalize construction of this interface | *----------------------------------------------------------------------*/ bool MOERTEL::Interface::Complete() { if (IsComplete()) { if (OutLevel()>0) std::cout << "MOERTEL: ***WRN*** MOERTEL::Interface::InterfaceComplete:\n" << "MOERTEL: ***WRN*** InterfaceComplete() was called before, do nothing\n" << "MOERTEL: ***WRN*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return true; } //------------------------------------------------------------------- // check for NULL entries in maps bool ok = true; for (int i=0; i<2; ++i) { std::map<int,Teuchos::RCP<MOERTEL::Node> >::const_iterator curr; for (curr=node_[i].begin(); curr!=node_[i].end(); ++curr) { if (curr->second == Teuchos::null) { std::cout << "***ERR*** MOERTEL::Interface::Complete:\n" << "***ERR*** Interface # " << Id_ << ":\n" << "***ERR*** found NULL entry in map of nodes\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; ok = false; } } } for (int i=0; i<2; ++i) { std::map<int,Teuchos::RCP<MOERTEL::Segment> >::const_iterator curr; for (curr=seg_[i].begin(); curr!=seg_[i].end(); ++curr) { if (curr->second == Teuchos::null) { std::cout << "***ERR*** MOERTEL::Interface::Complete:\n" << "***ERR*** Interface # " << Id_ << ":\n" << "***ERR*** found NULL entry in map of segments\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; ok = false; } } } int lok = ok; int gok = 1; gcomm_.MinAll(&lok,&gok,1); if (!gok) return false; //------------------------------------------------------------------- // check whether all nodes for segments are present // (take in account that node might be on different processor) // this test is expensive and does not scale. It is therefore only performed // when user requests a high output level #if 1 if (OutLevel()>9) { for (int proc=0; proc<gcomm_.NumProc(); ++proc) { for (int side=0; side<2; ++side) { // create length of list of all nodes adjacent to segments on proc int sendsize = 0; if (proc==gcomm_.MyPID()) { std::map<int,Teuchos::RCP<MOERTEL::Segment> >::const_iterator curr; for (curr=seg_[side].begin(); curr!=seg_[side].end(); ++curr) sendsize += curr->second->Nnode(); } gcomm_.Broadcast(&sendsize,1,proc); // create list of all nodes adjacent to segments on proc std::vector<int> ids(sendsize); if (proc==gcomm_.MyPID()) { std::map<int,Teuchos::RCP<MOERTEL::Segment> >::const_iterator curr; int counter=0; for (curr=seg_[side].begin(); curr!=seg_[side].end(); ++curr) { const int* segids = curr->second->NodeIds(); for (int i=0; i<curr->second->Nnode(); ++i) ids[counter++] = segids[i]; } } gcomm_.Broadcast(&ids[0],sendsize,proc); // check on all processors for nodes in ids std::vector<int> foundit(sendsize); std::vector<int> gfoundit(sendsize); for (int i=0; i<sendsize; ++i) { foundit[i] = 0; if (node_[side].find(ids[i]) != node_[side].end()) foundit[i] = 1; } gcomm_.MaxAll(&foundit[0],&gfoundit[0],sendsize); for (int i=0; i<sendsize; ++i) { if (gfoundit[i]!=1) { if (gcomm_.MyPID()==proc) std::cout << "***ERR*** MOERTEL::Interface::Complete:\n" << "***ERR*** cannot find segment's node # " << ids[i] << "\n" << "***ERR*** in map of all nodes on all procs\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; ids.clear(); foundit.clear(); gfoundit.clear(); gcomm_.Barrier(); return false; } } // tidy up ids.clear(); foundit.clear(); gfoundit.clear(); } // for (int size=0; side<2; ++side) } // for (int proc=0; proc<gcomm_.NumProc(); ++proc) } #endif //------------------------------------------------------------------- // find all procs that have business on this interface (own nodes/segments) // build a Epetra_comm that contains only those procs // this intra-communicator will be used to handle most stuff on this // interface so the interface will not block all other procs { #ifdef EPETRA_MPI std::vector<int> lin(gcomm_.NumProc()); std::vector<int> gin(gcomm_.NumProc()); for (int i=0; i<gcomm_.NumProc(); ++i) lin[i] = 0; // check ownership of any segments for (int i=0; i<2; ++i) if (seg_[i].size() != 0) { lin[gcomm_.MyPID()] = 1; break; } // check ownership of any nodes for (int i=0; i<2; ++i) if (node_[i].size() != 0) { lin[gcomm_.MyPID()] = 1; break; } gcomm_.MaxAll(&lin[0],&gin[0],gcomm_.NumProc()); lin.clear(); // typecast the Epetra_Comm to Epetra_MpiComm Epetra_MpiComm* epetrampicomm = dynamic_cast<Epetra_MpiComm*>(&gcomm_); if (!epetrampicomm) { std::stringstream oss; oss << "***ERR*** MOERTEL::Interface::Complete:\n" << "***ERR*** Interface " << Id() << ": Epetra_Comm is not an Epetra_MpiComm\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw ReportError(oss); } // split the communicator into participating and none-participating procs int color; int key = gcomm_.MyPID(); // I am taking part in the new comm if I have any ownership if (gin[gcomm_.MyPID()]) color = 0; // I am not taking part in the new comm else color = MPI_UNDEFINED; // tidy up gin.clear(); // create the local communicator MPI_Comm mpi_global_comm = epetrampicomm->GetMpiComm(); MPI_Comm* mpi_local_comm = new MPI_Comm(); MPI_Comm_split(mpi_global_comm,color,key,mpi_local_comm); // create the new Epetra_MpiComm if (*mpi_local_comm == MPI_COMM_NULL) lcomm_ = Teuchos::null; else lcomm_ = Teuchos::rcp(new Epetra_MpiComm(*mpi_local_comm)); // FIXME: who destroys the MPI_Comm inside? #if 0 // test this stuff on the mpi level int grank,lrank; MPI_Comm_rank(mpi_global_comm,&grank); if (*mpi_local_comm != MPI_COMM_NULL) MPI_Comm_rank(*mpi_local_comm,&lrank); else lrank = -1; for (int proc=0; proc<gcomm_.NumProc(); ++proc) { if (proc==gcomm_.MyPID()) std::cout << "using mpi comms: I am global rank " << grank << " and local rank " << lrank << std::endl; gcomm_.Barrier(); } // test this stuff on the epetra level if (lComm()) for (int proc=0; proc<lcomm_->NumProc(); ++proc) { if (proc==lcomm_->MyPID()) std::cout << "using epetra comms: I am global rank " << gcomm_.MyPID() << " and local rank " << lcomm_->MyPID() << std::endl; lcomm_->Barrier(); } gcomm_.Barrier(); #endif #else // the easy serial case Epetra_SerialComm* serialcomm = dynamic_cast<Epetra_SerialComm*>(&gcomm_); if (!serialcomm) { std::stringstream oss; oss << "***ERR*** MOERTEL::Interface::Complete:\n" << "***ERR*** Interface " << Id() << ": Epetra_Comm is not an Epetra_SerialComm\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw ReportError(oss); } lcomm_ = Teuchos::rcp(new Epetra_SerialComm(*serialcomm)); #endif // end of #ifdef PARALLEL } //------------------------------------------------------------------- // create a map of all nodes to there PID (process id) if (lComm()) for (int proc=0; proc<lcomm_->NumProc(); ++proc) { int lnnodes = 0; if (proc==lcomm_->MyPID()) lnnodes = node_[0].size() + node_[1].size(); lcomm_->Broadcast(&lnnodes,1,proc); std::vector<int> ids(lnnodes); if (proc==lcomm_->MyPID()) { std::map<int,Teuchos::RCP<MOERTEL::Node> >::const_iterator curr; int counter=0; for (int side=0; side<2; ++side) for (curr=node_[side].begin(); curr!=node_[side].end(); ++curr) ids[counter++] = curr->first; } lcomm_->Broadcast(&ids[0],lnnodes,proc); for (int i=0; i<lnnodes; ++i) nodePID_.insert(std::pair<int,int>(ids[i],proc)); ids.clear(); } //------------------------------------------------------------------- // create a map of all segments to there PID (process id) if (lComm()) for (int proc=0; proc<lcomm_->NumProc(); ++proc) { int lnsegs = 0; if (proc==lcomm_->MyPID()) lnsegs = seg_[0].size() + seg_[1].size(); lcomm_->Broadcast(&lnsegs,1,proc); std::vector<int> ids(lnsegs); if (proc==lcomm_->MyPID()) { std::map<int,Teuchos::RCP<MOERTEL::Segment> >::const_iterator curr; int counter=0; for (int side=0; side<2; ++side) for (curr=seg_[side].begin(); curr!=seg_[side].end(); ++curr) ids[counter++] = curr->first; } lcomm_->Broadcast(&ids[0],lnsegs,proc); for (int i=0; i<lnsegs; ++i) segPID_.insert(std::pair<int,int>(ids[i],proc)); ids.clear(); } //------------------------------------------------------------------- // set isComplete_ flag // we set it here already as we will be using some methods that require it // from now on isComplete_ = true; //------------------------------------------------------------------- // make the nodes know there adjacent segments // find max number of nodes to a segment if (lComm()) { int lmaxnnode = 0; int gmaxnnode = 0; for (int side=0; side<2; ++side) { std::map<int,Teuchos::RCP<MOERTEL::Segment> >::const_iterator scurr; for (scurr=seg_[side].begin(); scurr!=seg_[side].end(); ++scurr) if (lmaxnnode < scurr->second->Nnode()) lmaxnnode = scurr->second->Nnode(); } lcomm_->MaxAll(&lmaxnnode,&gmaxnnode,1); // loop all procs and broadcast their adjacency for (int proc=0; proc<lcomm_->NumProc(); ++proc) { // local number of segments int lnseg = 0; if (proc==lcomm_->MyPID()) lnseg = seg_[0].size() + seg_[1].size(); lcomm_->Broadcast(&lnseg,1,proc); // allocate vector to hold adjacency int offset = gmaxnnode+2; int size = lnseg*offset; std::vector<int> adj(size); // proc fills adjacency vector adj and broadcasts if (proc==lcomm_->MyPID()) { int count = 0; for (int side=0; side<2; ++side) { std::map<int,Teuchos::RCP<MOERTEL::Segment> >::const_iterator scurr; for (scurr=seg_[side].begin(); scurr!=seg_[side].end(); ++scurr) { Teuchos::RCP<MOERTEL::Segment> seg = scurr->second; adj[count] = seg->Id(); adj[count+1] = seg->Nnode(); const int* ids = seg->NodeIds(); for (int i=0; i<seg->Nnode(); ++i) adj[count+2+i] = ids[i]; count += offset; } } } lcomm_->Broadcast(&adj[0],size,proc); // all procs read adj and add segment to the nodes they own int count = 0; for (int i=0; i<lnseg; ++i) { int segid = adj[count]; int nnode = adj[count+1]; for (int j=0; j<nnode; ++j) { int nid = adj[count+2+j]; if (lcomm_->MyPID() == NodePID(nid)) { // I own this node, so set the segment segid in it Teuchos::RCP<MOERTEL::Node> node = GetNodeViewLocal(nid); if (node == Teuchos::null) { std::stringstream oss; oss << "***ERR*** MOERTEL::Interface::Complete:\n" << "***ERR*** cannot find node " << nid << "\n" << "***ERR*** in map of all nodes on this proc\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw ReportError(oss); } node->AddSegment(segid); } else continue; } count += offset; } adj.clear(); } // for (int proc=0; proc<lcomm_->NumProc(); ++proc) } // if (lComm()) //------------------------------------------------------------------- // build redundant segments and nodes if (lComm()) { int ok = 0; ok += RedundantSegments(0); ok += RedundantSegments(1); ok += RedundantNodes(0); ok += RedundantNodes(1); if (ok != 4) { std::stringstream oss; oss << "***ERR*** MOERTEL::Interface::Complete:\n" << "***ERR*** building of redundant information failed\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw ReportError(oss); } } //------------------------------------------------------------------- // make topology segments <-> nodes for each side if (lComm()) BuildNodeSegmentTopology(); //------------------------------------------------------------------- // delete distributed nodes and segments for (int i=0; i<2; ++i) { seg_[i].clear(); node_[i].clear(); } //------------------------------------------------------------------- // we are done // note that there might not be any functions on the interface yet // they still have to be set return ok; }
InterfaceT)::ProjectNodes_MastertoSlave_Orthogonal() { if (!IsComplete()) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_MastertoSlave_Orthogonal:\n" << "***ERR*** Complete() not called on interface " << Id() << "\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } if (lcomm_ == Teuchos::null) return true; int mside = MortarSide(); int sside = OtherSide(mside); // iterate over all master nodes and project those belonging to me std::map<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)>>::iterator mcurr; for (mcurr = rnode_[mside].begin(); mcurr != rnode_[mside].end(); ++mcurr) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> mnode = mcurr->second; if (NodePID(mnode->Id()) != lcomm_->getRank()) continue; const double* mx = mnode->XCoords(); double mindist = 1.0e+20; Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> closenode = Teuchos::null; // find a node on the slave side that is closest to me std::map<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)>>:: iterator scurr; for (scurr = rnode_[sside].begin(); scurr != rnode_[sside].end(); ++scurr) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> snode = scurr->second; const double* sx = snode->XCoords(); // build distance | snode->XCoords() - mnode->XCoords() | double dist = 0.0; for (int i = 0; i < 3; ++i) dist += (mx[i] - sx[i]) * (mx[i] - sx[i]); dist = sqrt(dist); if (dist < mindist) { mindist = dist; closenode = snode; } } if (closenode == Teuchos::null) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_MastertoSlave_Orthogonal:\n" << "***ERR*** Weired: for master node " << mnode->Id() << " no closest master node found\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } // get segments attached to closest node closenode int nseg = closenode->Nseg(); MoertelT::SEGMENT_TEMPLATE_CLASS(SegmentT)** segs = closenode->Segments(); // create a projection operator MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectorT) projector(IsOneDimensional(), OutLevel()); // loop these segments and find best projection double bestdist[2]; double bestgap = 0.0, gap; bestdist[0] = bestdist[1] = 1.0e+20; MoertelT::SEGMENT_TEMPLATE_CLASS(SegmentT)* bestseg = NULL; for (int i = 0; i < nseg; ++i) { // project the master node orthogonally on the slave segment double xi[2]; xi[0] = xi[1] = 0.0; projector.ProjectNodetoSegment_SegmentOrthogonal( *mnode, *(segs[i]), xi, gap); // check whether xi is better than previous projection if (IsOneDimensional()) { if (abs(xi[0]) < abs(bestdist[0])) { bestdist[0] = xi[0]; bestdist[1] = xi[1]; bestseg = segs[i]; bestgap = gap; } } else { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_MastertoSlave_Orthogonal:\n" << "***ERR*** not impl. for 3D\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } } // for (int i=0; i<nseg; ++i) // check whether this best projection is good bool ok = false; if (IsOneDimensional()) { if (abs(bestdist[0]) < 1.01) ok = true; } else { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_MastertoSlave_Orthogonal:\n" << "***ERR*** not impl. for 3D\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } if (ok) // the projection is good { // create a projected node and store it in mnode MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)* pnode = new MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)( *mnode, bestdist, bestseg); mnode->SetProjectedNode(pnode); mnode->SetGap(bestgap); } else // this mnode does not have a valid projection { if (OutLevel() > 6) std::cout << "MoertelT: ***WRN***: Node " << mnode->Id() << " does not have projection\n\n"; // mnode->SetProjectedNode(NULL); } } // for (mcurr=rnode_[mside].begin(); mcurr!=rnode_[mside].end(); ++mcurr) // loop all master nodes again and make projection redundant int bsize = 4 * rnode_[mside].size(); std::vector<double> bcast(bsize); for (int proc = 0; proc < lcomm_->getSize(); ++proc) { int blength = 0; if (proc == lcomm_->getRank()) { for (mcurr = rnode_[mside].begin(); mcurr != rnode_[mside].end(); ++mcurr) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> mnode = mcurr->second; if (proc != NodePID(mnode->Id())) continue; // cannot have a projection on a node i don't own Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)> pnode = mnode->GetProjectedNode(); if (pnode == Teuchos::null) continue; // this node does not have a projection const double* xi = pnode->Xi(); bcast[blength] = (double)pnode->Id(); ++blength; bcast[blength] = (double)pnode->Segment()->Id(); ++blength; bcast[blength] = xi[0]; ++blength; bcast[blength] = xi[1]; ++blength; bcast[blength] = pnode->Gap(); ++blength; } // for (mcurr=rnode_[mside].begin(); mcurr!=rnode_[mside].end(); // ++mcurr) if (blength > bsize) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_MastertoSlave_Orthogonal:\n" << "***ERR*** Overflow in communication buffer occured\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } } // if (proc==lComm()->MyPID()) Teuchos::broadcast<LO, int>(*lcomm_, proc, 1, &blength); Teuchos::broadcast<LO, double>(*lcomm_, proc, blength, &bcast[0]); if (proc != lcomm_->getRank()) { int i; for (i = 0; i < blength;) { int nid = (int)bcast[i]; ++i; int sid = (int)bcast[i]; ++i; double* xi = &bcast[i]; ++i; ++i; double gap = bcast[i]; ++i; Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> mnode = GetNodeView(nid); Teuchos::RCP<MoertelT::SEGMENT_TEMPLATE_CLASS(SegmentT)> seg = GetSegmentView(sid); if (mnode == Teuchos::null || seg == Teuchos::null) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_MastertoSlave_Orthogonal:\n" << "***ERR*** Cannot get view of node or segment\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)* pnode = new MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)( *mnode, xi, seg.get()); mnode->SetProjectedNode(pnode); mnode->SetGap(gap); } if (i != blength) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_MastertoSlave_Orthogonal:\n" << "***ERR*** Mismatch in dimension of recv buffer: " << i << " != " << blength << "\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } } // if (proc!=lComm()->MyPID()) } // for (int proc=0; proc<lComm()->NumProc(); ++proc) bcast.clear(); return true; }
InterfaceT)::ProjectNodes_SlavetoMaster_Orthogonal() { if (!IsComplete()) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_SlavetoMaster_Orthogonal:\n" << "***ERR*** Complete() not called on interface " << Id() << "\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } if (lcomm_ == Teuchos::null) return true; int mside = MortarSide(); int sside = OtherSide(mside); // iterate over all nodes of the slave side and project those belonging to me std::map<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)>>::iterator scurr; for (scurr = rnode_[sside].begin(); scurr != rnode_[sside].end(); ++scurr) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> snode = scurr->second; #if 0 std::cout << "now projecting\n " << *snode; #endif if (NodePID(snode->Id()) != lcomm_->getRank()) continue; const double* sx = snode->XCoords(); double mindist = 1.0e+20; Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> closenode = Teuchos::null; // find a node on the master side, that is closest to me std::map<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)>>:: iterator mcurr; for (mcurr = rnode_[mside].begin(); mcurr != rnode_[mside].end(); ++mcurr) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> mnode = mcurr->second; const double* mx = mnode->XCoords(); // build distance | mnode->XCoords() - snode->XCoords() | double dist = 0.0; for (int i = 0; i < 3; ++i) dist += (mx[i] - sx[i]) * (mx[i] - sx[i]); dist = sqrt(dist); if (dist <= mindist) { mindist = dist; closenode = mnode; } // std::cout << "snode " << snode->Id() << " mnode " << mnode->Id() << " // mindist " << mindist << " dist " << dist << std::endl; } if (closenode == Teuchos::null) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_SlavetoMaster_Orthogonal:\n" << "***ERR*** Weired: for slave node " << snode->Id() << " no closest master node found\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } // get segments attached to closest node cnode int nmseg = closenode->Nseg(); MoertelT::SEGMENT_TEMPLATE_CLASS(SegmentT)** msegs = closenode->Segments(); // create a projection-iterator MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectorT) projector(IsOneDimensional(), OutLevel()); // loop segments and find all projections onto them int nsseg = snode->Nseg(); MoertelT::SEGMENT_TEMPLATE_CLASS(SegmentT)** ssegs = snode->Segments(); for (int i = 0; i < nmseg; ++i) { // loop all segments that are adjacent to the slave node for (int j = 0; j < nsseg; ++j) { // project the slave node onto that master segment double xi[2]; xi[0] = xi[1] = 0.0; double gap; projector.ProjectNodetoSegment_Orthogonal_to_Slave( *snode, *(msegs[i]), xi, gap, *(ssegs[j])); // check whether this projection is good bool ok = false; if (IsOneDimensional()) { if (abs(xi[0]) < 1.01) ok = true; } else { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_SlavetoMaster_Orthogonal:\n" << "***ERR*** not impl. for 3D\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } if (ok) // the projection is good { // create a projected node and store it in snode MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)* pnode = new MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)( *snode, xi, msegs[i], ssegs[j]->Id()); snode->SetProjectedNode(pnode); snode->SetGap(gap); #if 0 std::cout << " snode id: " << pnode->Id() << " projects on mseg: " << msegs[i]->Id() << " orth to sseg " << ssegs[j]->Id() << std::endl; #endif } } // for (int j=0; j<nsseg; ++j) } // for (int i=0; i<nmseg; ++i) } // for (scurr=rnode_[sside].begin(); scurr!=rnode_[sside].end(); ++scurr) // loop all slave nodes again and make projections redundant if (lcomm_->getSize() > 1) { std::vector<double> bcast(10 * rnode_[sside].size()); for (int proc = 0; proc < lcomm_->getSize(); ++proc) { int blength = 0; if (proc == lcomm_->getRank()) { for (scurr = rnode_[sside].begin(); scurr != rnode_[sside].end(); ++scurr) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> snode = scurr->second; if (proc != NodePID(snode->Id())) continue; // cannot have a projection on a node i don't own int npnode = 0; Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)>* pnode = snode->GetProjectedNode(npnode); if (!pnode) continue; // no projection on this one bcast[blength] = (double)snode->Id(); ++blength; bcast[blength] = (double)npnode; ++blength; for (int j = 0; j < npnode; ++j) { bcast[blength] = (double)pnode[j]->Segment()->Id(); ++blength; const double* xi = pnode[j]->Xi(); bcast[blength] = xi[0]; ++blength; bcast[blength] = xi[1]; ++blength; bcast[blength] = pnode[j]->OrthoSegment(); ++blength; bcast[blength] = pnode[j]->Gap(); ++blength; } if ((int)bcast.size() < blength + 20) bcast.resize(bcast.size() + 40); } // for (mcurr=rnode_[mside].begin(); mcurr!=rnode_[mside].end(); // ++mcurr) if (blength >= (int)bcast.size()) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_SlavetoMaster_Orthogonal:\n" << "***ERR*** Overflow in communication buffer occured\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } } // if (proc==lComm()->MyPID()) Teuchos::broadcast<LO, int>(*lcomm_, proc, 1, &blength); if (proc != lcomm_->getRank()) bcast.resize(blength); Teuchos::broadcast<LO, double>(*lcomm_, proc, blength, &bcast[0]); if (proc != lcomm_->getRank()) { int i; for (i = 0; i < blength;) { int nid = (int)bcast[i]; ++i; Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> snode = GetNodeView(nid); int npnode = (int)bcast[i]; ++i; for (int j = 0; j < npnode; ++j) { int sid = (int)bcast[i]; ++i; double* xi = &bcast[i]; ++i; ++i; int orthseg = (int)bcast[i]; ++i; double gap = bcast[i]; ++i; Teuchos::RCP<MoertelT::SEGMENT_TEMPLATE_CLASS(SegmentT)> seg = GetSegmentView(sid); MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)* pnode = new MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)( *snode, xi, seg.get(), orthseg); snode->SetProjectedNode(pnode); snode->SetGap(gap); } } if (i != blength) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_SlavetoMaster_Orthogonal:\n" << "***ERR*** Mismatch in dimension of recv buffer: " << i << " != " << blength << "\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } } // if (proc!=lComm()->MyPID()) } // for (int proc=0; proc<lComm()->NumProc(); ++proc) bcast.clear(); } // if (lComm()->NumProc()>1) #if 1 // Postprocess the projections // The slave side of the interface might be larger then the master side // of the interface so not all slave nodes have a projection. // For those slave nodes without a projection attached to a slave segment // which overlaps with the master side, Lagrange multipliers have to be // introduced. This is done by checking all nodes without a projection // whether they are attached to some slave segment on which another node // HAS a projection. If this case is found, a pseudo ProjectedNode is // introduced for that node. for (scurr = rnode_[sside].begin(); scurr != rnode_[sside].end(); ++scurr) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> snode = scurr->second; // do only my own nodes // don't do anything on nodes that already have a projection if (snode->GetProjectedNode() != Teuchos::null) continue; // get segments adjacent to this node int nseg = snode->Nseg(); MoertelT::SEGMENT_TEMPLATE_CLASS(SegmentT)** segs = snode->Segments(); // loop segments and check for other nodes with projection bool foundit = false; for (int i = 0; i < nseg; ++i) { int nnode = segs[i]->Nnode(); MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)** nodes = segs[i]->Nodes(); for (int j = 0; j < nnode; ++j) if (nodes[j]->GetProjectedNode() != Teuchos::null) if (nodes[j]->GetProjectedNode()->Segment()) { foundit = true; break; } if (foundit) break; } if (foundit) { #if 0 std::cout << "Node without projection:\n" << *snode; std::cout << "...get's lagrange multipliers\n\n"; #endif MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)* pnode = new MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)( *snode, NULL, NULL); snode->SetProjectedNode(pnode); snode->SetGap(0.); } } #endif return true; }
InterfaceT)::ProjectNodes_SlavetoMaster_NormalField() { if (!IsComplete()) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_SlavetoMaster_NormalField:\n" << "***ERR*** Complete() not called on interface " << Id() << "\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } if (lcomm_ == Teuchos::null) return true; int mside = MortarSide(); int sside = OtherSide(mside); // iterate over all nodes of the slave side and project those belonging to me std::map<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)>>::iterator scurr; for (scurr = rnode_[sside].begin(); scurr != rnode_[sside].end(); ++scurr) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> snode = scurr->second; if (NodePID(snode->Id()) != lcomm_->getRank()) continue; const double* sx = snode->XCoords(); double mindist = 1.0e+20; Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> closenode = Teuchos::null; // find a node on the master side, that is closest to me std::map<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)>>:: iterator mcurr; for (mcurr = rnode_[mside].begin(); mcurr != rnode_[mside].end(); ++mcurr) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> mnode = mcurr->second; const double* mx = mnode->XCoords(); // build distance | mnode->XCoords() - snode->XCoords() | double dist = 0.0; for (int i = 0; i < 3; ++i) dist += (mx[i] - sx[i]) * (mx[i] - sx[i]); dist = sqrt(dist); if (dist <= mindist) { mindist = dist; closenode = mnode; } std::cout << "snode " << snode->Id() << " mnode " << mnode->Id() << " mindist " << mindist << " dist " << dist << std::endl; } if (closenode == Teuchos::null) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_SlavetoMaster_NormalField:\n" << "***ERR*** Weired: for slave node " << snode->Id() << " no closest master node found\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } #if 0 std::cout << "snode " << *snode; std::cout << "closenode " << *closenode; #endif // get segments attached to closest node cnode int nseg = closenode->Nseg(); MoertelT::SEGMENT_TEMPLATE_CLASS(SegmentT)** segs = closenode->Segments(); // create a projection-iterator MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectorT) projector(IsOneDimensional(), OutLevel()); // finding a good geometric projection is somehow // critical. We work with some tolerance here and pick the 'best' // out of all acceptable projections made // loop these segments and project onto them along snode's normal vector double bestdist[2]; double gap, bestgap = 0.0; const double tol = 0.2; bestdist[0] = bestdist[1] = 1.0e+20; MoertelT::SEGMENT_TEMPLATE_CLASS(SegmentT)* bestseg = NULL; for (int i = 0; i < nseg; ++i) { // project the slave node onto that master segment double xi[2]; xi[0] = xi[1] = 0.0; projector.ProjectNodetoSegment_NodalNormal(*snode, *(segs[i]), xi, gap); // check whether xi is better then previous projections if (IsOneDimensional()) // 2D case { if (abs(xi[0]) < abs(bestdist[0])) { bestdist[0] = xi[0]; bestdist[1] = xi[1]; bestseg = segs[i]; bestgap = gap; } } else // 3D case { double third = 1. / 3.; // it's 'inside' with some tolerance if (xi[0] <= 1. + tol && xi[1] <= abs(1. - xi[0]) + tol && xi[0] >= 0. - tol && xi[1] >= 0. - tol) { // it's better in both directions if (sqrt((xi[0] - third) * (xi[0] - third)) < sqrt((bestdist[0] - third) * (bestdist[0] - third)) && sqrt((xi[1] - third) * (xi[1] - third)) < sqrt((bestdist[1] - third) * (bestdist[1] - third))) { bestdist[0] = xi[0]; bestdist[1] = xi[1]; bestseg = segs[i]; bestgap = gap; } // it's better in one direction and 'in' in the other else if ( (sqrt((xi[0] - third) * (xi[0] - third)) < sqrt((bestdist[0] - third) * (bestdist[0] - third)) && xi[1] <= abs(1. - xi[0]) + tol && xi[1] >= 0. - tol) || (sqrt((xi[1] - third) * (xi[1] - third)) < sqrt((bestdist[1] - third) * (bestdist[1] - third)) && xi[0] <= 1. + tol && xi[0] >= 0. - tol)) { bestdist[0] = xi[0]; bestdist[1] = xi[1]; bestseg = segs[i]; bestgap = gap; } } } } // for (int i=0; i<nseg; ++i) // check whether the bestseg and bestdist are inside the segment // (with some tolerance of 20%) bool ok = false; if (IsOneDimensional()) { if (abs(bestdist[0]) < 1.2) ok = true; } else { if (bestdist[0] <= 1. + tol && bestdist[1] <= abs(1. - bestdist[0]) + tol && bestdist[0] >= 0. - tol && bestdist[1] >= 0. - tol) ok = true; } if (ok) // the projection is good { // create a projected node and store it in snode MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)* pnode = new MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)( *snode, bestdist, bestseg); snode->SetProjectedNode(pnode); snode->SetGap(bestgap); } else { if (OutLevel() > 6) std::cout << "MoertelT: ***WRN***: Projection s->m: Node " << snode->Id() << " does not have projection\n"; snode->SetProjectedNode(NULL); } } // for (scurr=rnode_[sside].begin(); scurr!=rnode_[sside].end(); ++scurr) lcomm_->barrier(); // loop all slave nodes again and make the projections redundant std::vector<double> bcast(4 * rnode_[sside].size()); for (int proc = 0; proc < lcomm_->getSize(); ++proc) { int blength = 0; if (proc == lcomm_->getRank()) { for (scurr = rnode_[sside].begin(); scurr != rnode_[sside].end(); ++scurr) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> snode = scurr->second; if (proc != NodePID(snode->Id())) continue; // I cannot have a projection on a node not owned by me Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)> pnode = snode->GetProjectedNode(); if (pnode == Teuchos::null) continue; // this node does not have a projection const double* xi = pnode->Xi(); bcast[blength] = (double)pnode->Id(); ++blength; if (pnode->Segment()) bcast[blength] = (double)pnode->Segment()->Id(); else bcast[blength] = -0.1; // indicating this node does not have // projection but lagrange multipliers ++blength; bcast[blength] = xi[0]; ++blength; bcast[blength] = xi[1]; ++blength; bcast[blength] = pnode->Gap(); ++blength; } if (blength > (int)(4 * rnode_[sside].size())) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_SlavetoMaster_NormalField:\n" << "***ERR*** Overflow in communication buffer occured\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } } Teuchos::broadcast<LO, int>(*lcomm_, proc, 1, &blength); Teuchos::broadcast<LO, double>(*lcomm_, proc, blength, &bcast[0]); if (proc != lcomm_->getRank()) { int i; for (i = 0; i < blength;) { int nid = (int)bcast[i]; ++i; double sid = bcast[i]; ++i; double* xi = &bcast[i]; ++i; ++i; double gap = bcast[i]; ++i; Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> snode = GetNodeView(nid); Teuchos::RCP<MoertelT::SEGMENT_TEMPLATE_CLASS(SegmentT)> seg = Teuchos::null; if (sid != -0.1) seg = GetSegmentView((int)sid); if (snode == Teuchos::null) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_SlavetoMaster_NormalField:" "\n" << "***ERR*** Cannot get view of node\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)* pnode = new MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)( *snode, xi, seg.get()); snode->SetProjectedNode(pnode); snode->SetGap(gap); } if (i != blength) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_SlavetoMaster_NormalField:\n" << "***ERR*** Mismatch in dimension of recv buffer\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } } } // for (int proc=0; proc<lComm()->NumProc(); ++proc) bcast.clear(); lcomm_->barrier(); #if 1 // Postprocess the projections // The slave side of the interface might be larger then the master side // of the interface so not all slave nodes have a projection. // For those slave nodes without a projection attached to a slave segment // which overlaps with the master side, Lagrange multipliers have to be // introduced. This is done by checking all nodes without a projection // whether they are attached to some slave segment on which another node // HAS a projection. If this case is found, a pseudo ProjectedNode is // introduced for that node. for (scurr = rnode_[sside].begin(); scurr != rnode_[sside].end(); ++scurr) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> snode = scurr->second; // don't do anything on nodes that already have a projection if (snode->GetProjectedNode() != Teuchos::null) continue; // get segments adjacent to this node int nseg = snode->Nseg(); MoertelT::SEGMENT_TEMPLATE_CLASS(SegmentT)** segs = snode->Segments(); // loop segments and check for other nodes with projection bool foundit = false; for (int i = 0; i < nseg; ++i) { int nnode = segs[i]->Nnode(); MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)** nodes = segs[i]->Nodes(); for (int j = 0; j < nnode; ++j) if (nodes[j]->GetProjectedNode() != Teuchos::null) if (nodes[j]->GetProjectedNode()->Segment()) { foundit = true; break; } if (foundit) break; } if (foundit) { #if 1 std::cout << "Node without projection:\n" << *snode; std::cout << "...get's lagrange multipliers\n\n"; #endif MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)* pnode = new MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)( *snode, NULL, NULL); snode->SetProjectedNode(pnode); snode->SetGap(0.); } } // for (scurr=rnode_[sside].begin(); scurr!=rnode_[sside].end(); ++scurr) #endif return true; }
InterfaceT)::ProjectNodes_MastertoSlave_NormalField() { if (!IsComplete()) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_MastertoSlave_NormalField:\n" << "***ERR*** Complete() not called on interface " << Id() << "\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } if (lcomm_ == Teuchos::null) return true; int mside = MortarSide(); int sside = OtherSide(mside); // iterate over all nodes of the master side and project those belonging to me std::map<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)>>::iterator mcurr; for (mcurr = rnode_[mside].begin(); mcurr != rnode_[mside].end(); ++mcurr) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> mnode = mcurr->second; if (NodePID(mnode->Id()) != lcomm_->getRank()) continue; const double* mx = mnode->XCoords(); double mindist = 1.0e+20; Teuchos::RCP<MoertelT::(NodeT)> closenode = Teuchos::null; // find a node on the slave side that is closest to me std::map<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)>>:: iterator scurr; for (scurr = rnode_[sside].begin(); scurr != rnode_[sside].end(); ++scurr) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> snode = scurr->second; const double* sx = snode->XCoords(); // build distance | snode->XCoords() - mnode->XCoords() | double dist = 0.0; for (int i = 0; i < 3; ++i) dist += (mx[i] - sx[i]) * (mx[i] - sx[i]); dist = sqrt(dist); if (dist < mindist) { mindist = dist; closenode = snode; } } if (closenode == Teuchos::null) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_MastertoSlave_NormalField:\n" << "***ERR*** Weired: for master node " << mnode->Id() << " no closest master node found\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } #if 0 std::cout << "snode " << *mnode; std::cout << "closenode " << *closenode; #endif // get segments attached to closest node closenode int nseg = closenode->Nseg(); MoertelT::SEGMENT_TEMPLATE_CLASS(SegmentT)** segs = closenode->Segments(); // create a projection operator MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectorT) projector(IsOneDimensional(), OutLevel()); // loop these segments and find best projection double bestdist[2]; const double tol = 0.2; double bestgap = 0.0, gap; bestdist[0] = bestdist[1] = 1.0e+20; MoertelT::SEGMENT_TEMPLATE_CLASS(SegmentT)* bestseg = NULL; for (int i = 0; i < nseg; ++i) { // project the master node on the slave segment along the segments // interpolated normal field double xi[2]; xi[0] = xi[1] = 0.0; projector.ProjectNodetoSegment_SegmentNormal(*mnode, *(segs[i]), xi, gap); // check whether xi is better then previous projections if (IsOneDimensional()) { if (abs(xi[0]) < abs(bestdist[0])) { bestdist[0] = xi[0]; bestdist[1] = xi[1]; bestseg = segs[i]; bestgap = gap; } } else { double third = 1. / 3.; // it's 'inside' with some tolerance if (xi[0] <= 1. + tol && xi[1] <= abs(1. - xi[0]) + tol && xi[0] >= 0. - tol && xi[1] >= 0. - tol) { // it's better in both directions if (sqrt((xi[0] - third) * (xi[0] - third)) < sqrt((bestdist[0] - third) * (bestdist[0] - third)) && sqrt((xi[1] - third) * (xi[1] - third)) < sqrt((bestdist[1] - third) * (bestdist[1] - third))) { bestdist[0] = xi[0]; bestdist[1] = xi[1]; bestseg = segs[i]; bestgap = gap; } // it's better in one direction and 'in' in the other else if ( (sqrt((xi[0] - third) * (xi[0] - third)) < sqrt((bestdist[0] - third) * (bestdist[0] - third)) && xi[1] <= abs(1. - xi[0]) + tol && xi[1] >= 0. - tol) || (sqrt((xi[1] - third) * (xi[1] - third)) < sqrt((bestdist[1] - third) * (bestdist[1] - third)) && xi[0] <= 1. + tol && xi[0] >= 0. - tol)) { bestdist[0] = xi[0]; bestdist[1] = xi[1]; bestseg = segs[i]; bestgap = gap; } } } } // for (int i=0; i<nseg; ++i) // check whether the bestseg/bestdist are inside that segment // (with some tolerance of 20%) bool ok = false; if (IsOneDimensional()) { if (abs(bestdist[0]) < 1.1) ok = true; } else { if (bestdist[0] <= 1. + tol && bestdist[1] <= abs(1. - bestdist[0]) + tol && bestdist[0] >= 0. - tol && bestdist[1] >= 0. - tol) ok = true; } if (ok) // the projection is good { // build the interpolated normal and overwrite the mnode normal with -n int nsnode = bestseg->Nnode(); MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)** snodes = bestseg->Nodes(); std::vector<double> val(nsnode); bestseg->EvaluateFunction(0, bestdist, &val[0], nsnode, NULL); double NN[3]; NN[0] = NN[1] = NN[2] = 0.0; for (int i = 0; i < nsnode; ++i) { const double* Normal = snodes[i]->Normal(); for (int j = 0; j < 3; ++j) NN[j] -= val[i] * Normal[j]; } val.clear(); mnode->SetN(NN); // create projected node and store it in mnode MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)* pnode = new MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)( *mnode, bestdist, bestseg); mnode->SetProjectedNode(pnode); mnode->SetGap(bestgap); } else // this mnode does not have a valid projection { if (OutLevel() > 6) std::cout << "MoertelT: ***WRN***: Projection m->s: Node " << mnode->Id() << " does not have projection\n"; mnode->SetProjectedNode(NULL); } } // for (scurr=rnode_[mside].begin(); scurr!=rnode_[mside].end(); ++scurr) // loop all master nodes again and make the projection and the new normal // redundant int bsize = 7 * rnode_[mside].size(); std::vector<double> bcast(bsize); for (int proc = 0; proc < lcomm_->getSize(); ++proc) { int blength = 0; if (proc == lcomm_->getRank()) { for (mcurr = rnode_[mside].begin(); mcurr != rnode_[mside].end(); ++mcurr) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> mnode = mcurr->second; if (proc != NodePID(mnode->Id())) continue; // cannot have a projection on a node i don't own Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)> pnode = mnode->GetProjectedNode(); if (pnode == Teuchos::null) continue; // this node does not have a projection const double* xi = pnode->Xi(); const double* Normal = mnode->Normal(); bcast[blength] = (double)pnode->Id(); ++blength; if (pnode->Segment()) bcast[blength] = (double)pnode->Segment()->Id(); else bcast[blength] = -0.1; ++blength; bcast[blength] = xi[0]; ++blength; bcast[blength] = xi[1]; ++blength; bcast[blength] = Normal[0]; ++blength; bcast[blength] = Normal[1]; ++blength; bcast[blength] = Normal[2]; ++blength; bcast[blength] = pnode->Gap(); ++blength; } if (blength > bsize) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_MastertoSlave_NormalField:\n" << "***ERR*** Overflow in communication buffer occured\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } } Teuchos::broadcast<LO, int>(*lcomm_, proc, 1, &blength); Teuchos::broadcast<LO, double>(*lcomm_, proc, blength, &bcast[0]); if (proc != lcomm_->getRank()) { int i; for (i = 0; i < blength;) { int nid = (int)bcast[i]; ++i; double sid = bcast[i]; ++i; double* xi = &bcast[i]; ++i; ++i; double* n = &bcast[i]; ++i; ++i; ++i; double gap = bcast[i]; ++i; Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(NodeT)> mnode = GetNodeView(nid); Teuchos::RCP<MoertelT::SEGMENT_TEMPLATE_CLASS(SegmentT)> seg = Teuchos::null; if (sid != -0.1) seg = GetSegmentView((int)sid); if (mnode == Teuchos::null) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_MastertoSlave_NormalField:" "\n" << "***ERR*** Cannot get view of node\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } mnode->SetN(n); MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)* pnode = new MoertelT::MOERTEL_TEMPLATE_CLASS(ProjectedNodeT)( *mnode, xi, seg.get()); mnode->SetProjectedNode(pnode); mnode->SetGap(gap); } if (i != blength) { std::stringstream oss; oss << "***ERR*** " "MoertelT::Interface::ProjectNodes_MastertoSlave_NormalField:\n" << "***ERR*** Mismatch in dimension of recv buffer\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw MoertelT::ReportError(oss); } } } // for (int proc=0; proc<lComm()->NumProc(); ++proc) bcast.clear(); return true; }
/*----------------------------------------------------------------------* | Create the saddle point problem (public) mwgee 07/05| | Note that this is collective for ALL procs | *----------------------------------------------------------------------*/ Epetra_CrsMatrix* MOERTEL::Manager::MakeSaddleProblem() { // time this process Epetra_Time time(Comm()); time.ResetStartTime(); // check whether all interfaces are complete and integrated std::map<int,Teuchos::RCP<MOERTEL::Interface> >::iterator curr; for (curr=interface_.begin(); curr != interface_.end(); ++curr) { if (curr->second->IsComplete() == false) { cout << "***ERR*** MOERTEL::Manager::MakeSaddleProblem:\n" << "***ERR*** interface " << curr->second->Id() << " is not Complete()\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } if (curr->second->IsIntegrated() == false) { cout << "***ERR*** MOERTEL::Manager::MakeSaddleProblem:\n" << "***ERR*** interface " << curr->second->Id() << " is not integrated yet\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } } // check whether we have a problemmap_ if (problemmap_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSaddleProblem:\n" << "***ERR*** No problemmap_ set\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // check whether we have a constraintsmap_ if (constraintsmap_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSaddleProblem:\n" << "***ERR*** onstraintsmap is NULL\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // check for saddlemap_ if (saddlemap_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSaddleProblem:\n" << "***ERR*** saddlemap_==NULL\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // check for inputmatrix if (inputmatrix_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSaddleProblem:\n" << "***ERR*** No inputmatrix set\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // check whether we have M and D matrices if (D_==Teuchos::null || M_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSaddleProblem:\n" << "***ERR*** Matrix M or D is NULL\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // create a matrix for the saddle problem and fill it saddlematrix_ = Teuchos::rcp(new Epetra_CrsMatrix(Copy,*saddlemap_,90)); // add values from inputmatrix MOERTEL::MatrixMatrixAdd(*inputmatrix_,false,1.0,*saddlematrix_,0.0); // add values from D_ MOERTEL::MatrixMatrixAdd(*D_,false,1.0,*saddlematrix_,1.0); MOERTEL::MatrixMatrixAdd(*D_,true,1.0,*saddlematrix_,1.0); // add values from M_ MOERTEL::MatrixMatrixAdd(*M_,false,1.0,*saddlematrix_,1.0); MOERTEL::MatrixMatrixAdd(*M_,true,1.0,*saddlematrix_,1.0); saddlematrix_->FillComplete(); saddlematrix_->OptimizeStorage(); // time this process double t = time.ElapsedTime(); if (OutLevel()>5 && Comm().MyPID()==0) cout << "MOERTEL (Proc 0): Construct saddle system in " << t << " sec\n"; return saddlematrix_.get(); }
/*----------------------------------------------------------------------* | Create the spd system (public) mwgee 12/05| | Note that this is collective for ALL procs | *----------------------------------------------------------------------*/ Epetra_CrsMatrix* MOERTEL::Manager::MakeSPDProblem() { // time this process Epetra_Time time(Comm()); time.ResetStartTime(); // check whether all interfaces are complete and integrated std::map<int,Teuchos::RCP<MOERTEL::Interface> >::iterator curr; for (curr=interface_.begin(); curr != interface_.end(); ++curr) { if (curr->second->IsComplete() == false) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** interface " << curr->second->Id() << " is not Complete()\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } if (curr->second->IsIntegrated() == false) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** interface " << curr->second->Id() << " is not integrated yet\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } } // check whether we have a problemmap_ if (problemmap_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** No problemmap_ set\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // check whether we have a constraintsmap_ if (constraintsmap_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** onstraintsmap is NULL\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // check for saddlemap_ if (saddlemap_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** saddlemap_==NULL\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // check for inputmatrix if (inputmatrix_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** No inputmatrix set\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // check whether we have M and D matrices if (D_==Teuchos::null || M_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** Matrix M or D is NULL\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // we need a map from lagrange multiplier dofs to primal dofs on the same node std::vector<MOERTEL::Node*> nodes(0); std::map<int,int> lm_to_dof; for (curr=interface_.begin(); curr!=interface_.end(); ++curr) { Teuchos::RCP<MOERTEL::Interface> inter = curr->second; inter->GetNodeView(nodes); for (int i=0; i<(int)nodes.size(); ++i) { if (!nodes[i]->Nlmdof()) continue; const int* dof = nodes[i]->Dof(); const int* lmdof = nodes[i]->LMDof(); for (int j=0; j<nodes[i]->Nlmdof(); ++j) { //cout << "j " << j << " maps lmdof " << lmdof[j] << " to dof " << dof[j] << endl; lm_to_dof[lmdof[j]] = dof[j]; } } } lm_to_dof_ = Teuchos::rcp(new std::map<int,int>(lm_to_dof)); // this is a very useful map for the moertel_ml_preconditioner /* _ _ | | | Arr Arn Mr | S = | | | Anr Ann D | | | MrT D 0 | |_ _ | _ _ | | | Arr Arn | A = | | | Anr Ann | |_ _| 1) Ann is square and we need it's Range/DomainMap annmap _ _ WT = |_ 0 Dinv _| 2) Build WT (has rowmap/rangemap annmap and domainmap problemmap_) _ _ | | | Mr | B = | | | D | |_ _| 3) Build B (has rowmap/rangemap problemmap_ and domainmap annmap) 4) Build I, the identity matrix with maps problemmap_,problemmap_); */ int err=0; //-------------------------------------------------------------------------- // 1) create the rangemap of Ann std::vector<int> myanngids(problemmap_->NumMyElements()); int count=0; std::map<int,int>::iterator intintcurr; for (intintcurr=lm_to_dof.begin(); intintcurr!=lm_to_dof.end(); ++intintcurr) { if (problemmap_->MyGID(intintcurr->second)==false) continue; if ((int)myanngids.size()<=count) myanngids.resize(myanngids.size()+50); myanngids[count] = intintcurr->second; ++count; } myanngids.resize(count); int numglobalelements; Comm().SumAll(&count,&numglobalelements,1); Epetra_Map* annmap = new Epetra_Map(numglobalelements,count,&myanngids[0],0,Comm()); annmap_ = Teuchos::rcp(annmap); myanngids.clear(); //-------------------------------------------------------------------------- // 2) create WT Epetra_CrsMatrix* WT = new Epetra_CrsMatrix(Copy,*annmap,1,false); for (intintcurr=lm_to_dof.begin(); intintcurr!=lm_to_dof.end(); ++intintcurr) { int lmdof = intintcurr->first; int dof = intintcurr->second; if (D_->MyGRID(lmdof)==false) continue; int lmlrid = D_->LRID(lmdof); int numentries; int* indices; double* values; err = D_->ExtractMyRowView(lmlrid,numentries,values,indices); if (err) cout << "D_->ExtractMyRowView returned err=" << err << endl; bool foundit = false; for (int j=0; j<numentries; ++j) { int gcid = D_->GCID(indices[j]); if (gcid<0) cout << "Cannot find gcid for indices[j]\n"; //cout << "Proc " << Comm().MyPID() << " lmdof " << lmdof << " dof " << dof << " gcid " << gcid << " val " << values[j] << endl; if (gcid==dof) { double val = 1./values[j]; err = WT->InsertGlobalValues(dof,1,&val,&dof); if (err<0) cout << "WT->InsertGlobalValues returned err=" << err << endl; foundit = true; break; } } if (!foundit) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** Cannot compute inverse of D_\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; cout << "lmdof " << lmdof << " dof " << dof << endl; return NULL; } } WT->FillComplete(*problemmap_,*annmap); WT_ = Teuchos::rcp(WT); //-------------------------------------------------------------------------- // 3) create B // create a temporary matrix with rowmap of the Ann block Epetra_CrsMatrix* tmp = new Epetra_CrsMatrix(Copy,*annmap,120); std::vector<int> newindices(100); for (intintcurr=lm_to_dof.begin(); intintcurr!=lm_to_dof.end(); ++intintcurr) { int lmdof = intintcurr->first; int dof = intintcurr->second; if (D_->MyGRID(lmdof)==false) continue; int lmlrid = D_->LRID(lmdof); if (lmlrid<0) cout << "Cannot find lmlrid for lmdof\n"; int numentries; int* indices; double* values; // extract and add values from D err = D_->ExtractMyRowView(lmlrid,numentries,values,indices); if (err) cout << "D_->ExtractMyRowView returned err=" << err << endl; if (numentries>(int)newindices.size()) newindices.resize(numentries); for (int j=0; j<numentries; ++j) { newindices[j] = D_->GCID(indices[j]); if (newindices[j]<0) cout << "Cannot find gcid for indices[j]\n"; } //cout << "Inserting from D in row " << dof << " cols/val "; //for (int j=0; j<numentries; ++j) cout << newindices[j] << "/" << values[j] << " "; //cout << endl; err = tmp->InsertGlobalValues(dof,numentries,values,&newindices[0]); if (err) cout << "tmp->InsertGlobalValues returned err=" << err << endl; // extract and add values from M err = M_->ExtractMyRowView(lmlrid,numentries,values,indices); if (err) cout << "M_->ExtractMyRowView returned err=" << err << endl; if (numentries>(int)newindices.size()) newindices.resize(numentries); for (int j=0; j<numentries; ++j) { newindices[j] = M_->GCID(indices[j]); if (newindices[j]<0) cout << "Cannot find gcid for indices[j]\n"; } //cout << "Inserting from M in row " << dof << " cols/val "; //for (int j=0; j<numentries; ++j) cout << newindices[j] << "/" << values[j] << " "; //cout << endl; err = tmp->InsertGlobalValues(dof,numentries,values,&newindices[0]); if (err) cout << "tmp->InsertGlobalValues returned err=" << err << endl; } tmp->FillComplete(*(problemmap_.get()),*annmap); // B is transposed of tmp EpetraExt::RowMatrix_Transpose* trans = new EpetraExt::RowMatrix_Transpose(false); Epetra_CrsMatrix* B = &(dynamic_cast<Epetra_CrsMatrix&>(((*trans)(const_cast<Epetra_CrsMatrix&>(*tmp))))); delete tmp; tmp = NULL; B_ = Teuchos::rcp(new Epetra_CrsMatrix(*B)); newindices.clear(); //-------------------------------------------------------------------------- // 4) create I Epetra_CrsMatrix* I = MOERTEL::PaddedMatrix(*problemmap_,1.0,1); I->FillComplete(*problemmap_,*problemmap_); I_ = Teuchos::rcp(I); //-------------------------------------------------------------------------- // 5) Build BWT = B * WT Epetra_CrsMatrix* BWT = MOERTEL::MatMatMult(*B,false,*WT,false,OutLevel()); //-------------------------------------------------------------------------- // 6) create CBT = (I-BWT)*A*W*BT Epetra_CrsMatrix* spdrhs = new Epetra_CrsMatrix(Copy,*problemmap_,10,false); MOERTEL::MatrixMatrixAdd(*BWT,false,-1.0,*spdrhs,0.0); MOERTEL::MatrixMatrixAdd(*I,false,1.0,*spdrhs,1.0); spdrhs->FillComplete(); spdrhs->OptimizeStorage(); spdrhs_ = Teuchos::rcp(spdrhs); Epetra_CrsMatrix* tmp1 = MOERTEL::MatMatMult(*inputmatrix_,false,*WT,true,OutLevel()); Epetra_CrsMatrix* tmp2 = MOERTEL::MatMatMult(*tmp1,false,*B,true,OutLevel()); delete tmp1; Epetra_CrsMatrix* CBT = MOERTEL::MatMatMult(*spdrhs,false,*tmp2,false,OutLevel()); delete tmp2; //-------------------------------------------------------------------------- // 7) create spdmatrix_ = A - CBT - (CBT)^T spdmatrix_ = Teuchos::rcp(new Epetra_CrsMatrix(Copy,*problemmap_,10,false)); MOERTEL::MatrixMatrixAdd(*inputmatrix_,false,1.0,*spdmatrix_,0.0); MOERTEL::MatrixMatrixAdd(*CBT,false,-1.0,*spdmatrix_,1.0); MOERTEL::MatrixMatrixAdd(*CBT,true,-1.0,*spdmatrix_,1.0); delete CBT; CBT = NULL; spdmatrix_->FillComplete(); spdmatrix_->OptimizeStorage(); //-------------------------------------------------------------------------- // tidy up lm_to_dof.clear(); delete trans; B = NULL; // time this process double t = time.ElapsedTime(); if (OutLevel()>5 && Comm().MyPID()==0) cout << "MOERTEL (Proc 0): Construct spd system in " << t << " sec\n"; return spdmatrix_.get(); }
/*----------------------------------------------------------------------* | unpack all data from a vector into this class mwgee 10/05| *----------------------------------------------------------------------*/ bool MOERTEL::Segment_BiLinearTri::UnPack(int* pack) { // note: first there has to be the size and second there has to be the type int count = 0; int size = pack[count++]; stype_ = (MOERTEL::Segment::SegmentType)pack[count++]; Id_ = pack[count++]; nodeId_.resize(pack[count++]); for (int i=0; i<(int)nodeId_.size(); ++i) nodeId_[i] = pack[count++]; int nfunc = pack[count++]; for (int i=0; i<nfunc; ++i) { int id = pack[count++]; int type = pack[count++]; MOERTEL::Function* func = MOERTEL::AllocateFunction((MOERTEL::Function::FunctionType)type,OutLevel()); Teuchos::RCP<MOERTEL::Function> rcptrfunc = Teuchos::rcp(func); functions_.insert(std::pair<int,Teuchos::RCP<MOERTEL::Function> >(id,rcptrfunc)); } if (count != size) { std::stringstream oss; oss << "***ERR*** MOERTEL::Segment_BiLinearTri::UnPack:\n" << "***ERR*** mismatch in packing size\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw ReportError(oss); } return true; }
/*----------------------------------------------------------------------* | perform clipping algorithm (private) mwgee 10/05| *----------------------------------------------------------------------*/ bool MOERTEL::Overlap::Clipelements() { if (!havemxi_ || !havesxi_ || !havelines_ || !havesxim_ || !havelinem_) { std::cout << "***ERR*** MOERTEL::Overlap::Clipelements:\n" << "***ERR*** initialization of Overlap class missing\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; exit(EXIT_FAILURE); } const int nmnode = mseg_.Nnode(); const int nsnode = sseg_.Nnode(); // put all mseg nodes in polygon for (int i=0; i<nmnode; ++i) AddPointtoPolygon(100+i,&mline_[i][0]); //for (int i=0; i<nmnode; ++i) // std::cout << "mline_[" << i << "][0] " << mline_[i][0] << " mline_[" << i << "][1] " << mline_[i][1] // << " mline_[" << i << "][2] " << mline_[i][2] << " mline_[" << i << "][3] " << mline_[i][3] << endl; //=========================================================================== // loop edges of slave segment and clip the master edges // add all intersection points that can be found for (int clipedge=0; clipedge<nsnode; ++clipedge) { // point on that clip edge (dim 2) double* PE = &sline_[clipedge][0]; // the outward normal to the clip edge (dim 2) double* N = &sn_[clipedge][0]; // loop edges of master segment and clip them against the clip edge // add all intersections for (int medge=0; medge<nmnode; ++medge) { bool ok; // start point of medge with id 100+medge double* P0 = &mline_[medge][0]; //int id0 = 100+medge; // end point of medge has id 100+(medge+1 < 3) double* P1 = &mline_[medge][2]; // find intersection between medge and clipedge // if there's an intersection, put it in the polygon // the id is 10*clipedge+medge double xi[2]; ok = Clip_Intersect(N,PE,P0,P1,xi); if (ok) { //std::cout << "OVERLAP Clipelements: adding intersection point " << 10*clipedge+medge << endl; AddPointtoPolygon(10*clipedge+medge,xi); } } // for (int medge=0; medge<3; ++medge) } // for (int clipedge=0; clipedge<3; ++clipedge) //=========================================================================== // loop all clipedges and clip all points incl intersections // that are in the polygon // throw away all points that are not inside { int np = SizePointPolygon(); std::vector<Teuchos::RCP<MOERTEL::Point> > point; PointView(point); int p; for (p=0; p<np; ++p) { bool ok = true; //std::cout << "OVERLAP Clipelements: now clipping point " << point[p]->Id() << endl; const double* P = point[p]->Xi(); for (int clipedge=0; clipedge<nsnode; ++clipedge) { // point on that clip edge (dim 2) double* PE = &sline_[clipedge][0]; // the outward normal to the clip edge (dim 2) double* N = &sn_[clipedge][0]; // clip point P against this edge ok = Clip_TestPoint(N,PE,P,1.0e-10); // leave point in if (ok) continue; // remove this point else { ok = false; break; } } // for (int clipedge=0; clipedge<3; ++clipedge) // if not found inside, remove point if (!ok) { //std::cout << "OVERLAP Clipelements: removing point " << point[p]->Id() << endl; RemovePointfromPolygon(point[p]->Id(),P); } } // for (int p=0; p<np; ++p) point.clear(); } //=========================================================================== // loop all slave nodes and clip against master element. // put those slave nodes that are inside master element into polygon // Note that this works in master segment coords and the node that is put in // is put in with slave segment coords as the polygon is completely in // slave segment coords // master point have ids 100,101,102 // slave points have ids 1000,1001,1002 // edge intersections have ids sedge*10+medge // try only to put slave points in if the polygon is not empty int np = SizePointPolygon(); if (np) { for (int i=0; i<nsnode; ++i) { bool ok = true; // get the slave point in master coords double* P = sxim_[i]; // loop master clip edges for (int clipedge=0; clipedge<nmnode; ++clipedge) { // point on master clipedge double* PE = &mlinem_[clipedge][0]; // the outward normal to the clip edge (dim 2) double* N = &mn_[clipedge][0]; // clip point P against this edge ok = Clip_TestPoint(N,PE,P,1.0e-5); // put point in if (ok) continue; else { ok = false; break; } } // for (int clipedge=0; clipedge<3; ++clipedge) // don't put point in if (!ok) continue; else { //std::cout << "OVERLAP Clipelements: inserting slave point " << 1000+i << " xi=" // << sxi_[i][0] << "/" << sxi_[i][1] << endl; AddPointtoPolygon(1000+i,sxi_[i]); } } // for (int i=0; i<3; ++i) } //=========================================================================== //=========================================================================== //=========================================================================== #if 0 // make printout of the polygon so far { int np = SizePointPolygon(); std::vector<Teuchos::RCP<MOERTEL::Point> > point; PointView(point); for (int p=0; p<np; ++p) { std::cout << "OVERLAP Clipelements: point " << setw(3) << point[p]->Id() << " xi " << point[p]->Xi()[0] << "/" << point[p]->Xi()[1] << endl; } point.clear(); } #endif //=========================================================================== // count how many corner nodes of mseg are in and how many // intersections there are np = SizePointPolygon(); //=========================================================================== // if there are no points, there is still the chance that all slave points // are inside the master element if (!np) { for (int i=0; i<3; ++i) { bool ok = true; // get the slave point in master coords double* P = sxim_[i]; for (int clipedge=0; clipedge<3; ++clipedge) { // point on master clipedge double* PE = &mlinem_[clipedge][0]; // the outward normal to the clip edge (dim 2) double* N = &mn_[clipedge][0]; // clip point P against this edge ok = Clip_TestPoint(N,PE,P,1.0e-5); // put point in if (ok) continue; else { ok = false; break; } } // for (int clipedge=0; clipedge<3; ++clipedge) // don't put point in if (!ok) continue; else { //std::cout << "OVERLAP Clipelements: inserting slave point " << 1000+i << " xi=" // << sxi_[i][0] << "/" << sxi_[i][1] << endl; AddPointtoPolygon(1000+i,sxi_[i]); } } // for (int i=0; i<3; ++i) //========================================================================= // check again how many points there are inside np = SizePointPolygon(); if (np>2); // all slave points are in master segment, just continue else if (np && np <= 2) { if (inter_.OutLevel()>8) std::cout << "MOERTEL: ***WRN*** MOERTEL::Overlap::Clipelements:\n" << "MOERTEL: ***WRN*** " << np << " slave nodes seem to be in master segment but no overlap detected\n" << "MOERTEL: ***WRN*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } else // with none or less then 3 points in we assume no overlap return false; } //=========================================================================== // maybe its a good idea to collapse intersections with nodes (that are in) // that are really close to each other if (p_.size()>3) CollapsePoints(p_,1.0e-5); #if 0 // make printout of the polygon so far { std::cout << "--------------------------------------------\n"; int np = SizePointPolygon(); std::vector<Teuchos::RCP<MOERTEL::Point> > point; PointView(point); for (int p=0; p<np; ++p) { std::cout << "OVERLAP Clipelements: point " << setw(3) << point[p]->Id() << " xi " << point[p]->Xi()[0] << "/" << point[p]->Xi()[1] << endl; } point.clear(); } #endif //=========================================================================== // check again np = SizePointPolygon(); if (np && np<3) { if (OutLevel()>8) std::cout << "MOERTEL: ***WRN*** MOERTEL::Overlap::Clipelements:\n" << "MOERTEL: ***WRN*** " << np << " nodes in polygon but could not detect overlap\n" << "MOERTEL: ***WRN*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } else if (!np) return false; //=========================================================================== // finish the polygon ConvexHull(p_); return true; }
/*----------------------------------------------------------------------* | mwgee 07/05| | modded by gah 07/2010 | *----------------------------------------------------------------------*/ bool MOERTEL::Projector::ProjectNodetoSegment_NodalNormal(MOERTEL::Node& node, MOERTEL::Segment& seg, double xi[], double &gap) { bool ok = true; // 2D version of the problem if (IsTwoDimensional()) { #if 0 std::cout << "----- Projector: Node " << node.Id() << " Segment " << seg.Id() << std::endl; std::cout << "Segment\n" << seg; MOERTEL::Node** nodes = seg.Nodes(); std::cout << *nodes[0]; std::cout << *nodes[1]; #endif // we do a newton iteration for the projection coordinates xi // set starting value to the middle of the segment double eta = 0.0; int i = 0; double F=0.0,dF=0.0,deta=0.0; for (i=0; i<10; ++i) { F = evaluate_F_2D_NodalNormal(node,seg,eta,gap); if (abs(F) < 1.0e-10) break; dF = evaluate_gradF_2D_NodalNormal(node,seg,eta); deta = (-F)/dF; eta += deta; } if (abs(F)>1.0e-9) { ok = false; if (OutLevel()>3) std::cout << "MOERTEL: ***WRN*** MOERTEL::Projector::ProjectNodetoSegment_NodalNormal:\n" << "MOERTEL: ***WRN*** Newton iteration failed to converge\n" << "MOERTEL: ***WRN*** #iterations = " << i << std::endl << "MOERTEL: ***WRN*** F(eta) = " << F << " gradF(eta) = " << dF << " eta = " << eta << " delta(eta) = " << deta << "\n" << "MOERTEL: ***WRN*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; } #if 0 std::cout << "#iterations = " << i << " F = " << F << " eta = " << eta << std::endl; #endif xi[0] = eta; return ok; } // 3D version of the problem else { #if 0 std::cout << "----- Projector: Node " << node.Id() << " Segment " << seg.Id() << std::endl; std::cout << "Segment " << seg; MOERTEL::Node** nodes = seg.Nodes(); std::cout << *nodes[0]; std::cout << *nodes[1]; std::cout << *nodes[2]; #endif // we do a newton iteration for the projection coordinates xi // set starting value to the middle of the segment double eta[2]; if (seg.Nnode()==3) eta[0] = eta[1] = 1./3.; else eta[0] = eta[1] = 0.; double alpha = 0.0; int i=0; double F[3], dF[3][3], deta[3]; double eps; for (i=0; i<30; ++i) { evaluate_FgradF_3D_NodalNormal(F,dF,node,seg,eta,alpha,gap); eps = MOERTEL::dot(F,F,3); if (eps < 1.0e-10) break; // std::cout << eps << std::endl; MOERTEL::solve33(dF,deta,F); eta[0] -= deta[0]; eta[1] -= deta[1]; alpha -= deta[2]; } if (eps>1.0e-10) { ok = false; if (OutLevel()>3) std::cout << "MOERTEL: ***WRN*** MOERTEL::Projector::ProjectNodetoSegment_NodalNormal:\n" << "MOERTEL: ***WRN*** 3D Newton iteration failed to converge\n" << "MOERTEL: ***WRN*** #iterations = " << i << std::endl << "MOERTEL: ***WRN*** eps = " << eps << " eta[3] = " << eta[0] << "/" << eta[1] << "/" << alpha << "\n" << "MOERTEL: ***WRN*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; } #if 0 if (i>10) std::cout << "#iterations = " << i << " eps = " << eps << " eta = " << eta[0] << "/" << eta[1] << std::endl; #endif xi[0] = eta[0]; xi[1] = eta[1]; } // return ok; }
bool MoertelT::MOERTEL_TEMPLATE_CLASS_1A(OverlapT, IFace)::ConvexHull( std::map<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(PointT)>>& p) { // # points int np = p.size(); // sort points by xi[0] coordinate std::vector<Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(PointT)>> points; PointView(p, points); std::vector<double> dlistx(np); std::vector<double> dlisty(np); std::vector<int> list2(np); for (int i = 0; i < np; ++i) { dlistx[i] = points[i]->Xi()[0]; list2[i] = i; } MOERTEL::sort(&dlistx[0], np, &list2[0]); // get assoc. y-list for (int i = 0; i < np; ++i) dlisty[i] = points[list2[i]]->Xi()[1]; // sort list again for those points where x is identical for (int i = 0; i < np - 1; ++i) if (dlistx[i] == dlistx[i + 1]) if (dlisty[i] > dlisty[i + 1]) { MOERTEL::swap(dlistx[i], dlistx[i + 1]); MOERTEL::swap(dlisty[i], dlisty[i + 1]); MOERTEL::swap(list2[i], list2[i + 1]); } // create a new polygon and put points in in sorted order std::map<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(PointT)>> newp; for (int i = 0; i < np; ++i) { Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(PointT)> tmp = Teuchos::rcp(new MoertelT::MOERTEL_TEMPLATE_CLASS(PointT)( i, points[list2[i]]->Xi(), OutLevel())); newp.insert( std::pair<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(PointT)>>( i, tmp)); } // delete the old one p.clear(); // copy new one over CopyPointPolygon(newp, p); // destroy the new one newp.clear(); points.clear(); std::map<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(PointT)>>:: iterator pcurr; #if 0 // printout the polygon std::cout << "Input polygon:\n"; for (pcurr=p.begin(); pcurr != p.end(); ++pcurr) if (pcurr->second != Teuchos::null) std::cout << *(pcurr->second); #endif //=========================================================================== // build the upper hull std::map<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(PointT)>> upper; PointView(p, points); // put in the first 2 points AddPointtoPolygon(upper, 0, points[0]->Xi()); // std::cout << *points[0]; AddPointtoPolygon(upper, 1, points[1]->Xi()); // std::cout << *points[1]; //--------------------------------------------------------------------------- for (int i = 2; i < np; ++i) { // add point[i] to upper AddPointtoPolygon(upper, i, points[i]->Xi()); // std::cout << *points[i]; // find whether we still have a convex hull while (upper.size() > 2 && !MakeRightTurnUpper(i, upper)) RemovePointBefore(i, upper); #if 0 // printout the current upper hull std::map<int,MoertelT::MOERTEL_TEMPLATE_CLASS(PointT)*>::iterator pcurr; for (pcurr=upper.begin(); pcurr != upper.end(); ++pcurr) if (pcurr->second) std::cout << *(pcurr->second); #endif } // for (int i=2; i<np; ++i) //=========================================================================== // build the lower hull std::map<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(PointT)>> lower; // put in the first 2 points AddPointtoPolygon( lower, np - 1, points[np - 1]->Xi()); // std::cout << *points[np-1]; AddPointtoPolygon( lower, np - 2, points[np - 2]->Xi()); // std::cout << *points[np-2]; //--------------------------------------------------------------------------- for (int i = np - 3; i >= 0; --i) { // add point[i] to lower AddPointtoPolygon(lower, i, points[i]->Xi()); // std::cout << *points[i]; // find whether we still have a convex hull while (lower.size() > 2 && !MakeRightTurnLower(i, lower)) RemovePointAfter(i, lower); #if 0 // printout the current lower hull std::map<int,Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(PointT) > >::iterator pcurr; for (pcurr=lower.begin(); pcurr != lower.end(); ++pcurr) if (pcurr->second) std::cout << *(pcurr->second); #endif } // for (int i=np-3; i>=0; --i) //=========================================================================== points.clear(); //=========================================================================== // join upper and lower hull // note not to put in duplicate start and end point std::map<int, Teuchos::RCP<MoertelT::MOERTEL_TEMPLATE_CLASS(PointT)>> finalp; // put upper hull in int i = 0; for (pcurr = upper.begin(); pcurr != upper.end(); ++pcurr) { // std::cout << *(pcurr->second); AddPointtoPolygon(finalp, i, pcurr->second->Xi()); ++i; } // put lower hull in, skip first and last point pcurr = lower.end(); --pcurr; --pcurr; for (; pcurr != lower.begin(); --pcurr) { // std::cout << *(pcurr->second); AddPointtoPolygon(finalp, i, pcurr->second->Xi()); ++i; } #if 0 // printout the polygon std::cout << "--------------------------------------------\n"; std::cout << "Final polygon:\n"; for (pcurr=finalp.begin(); pcurr != finalp.end(); ++pcurr) if (pcurr->second != Teuchos::null) std::cout << *(pcurr->second); #endif // compare size of convex hull polygon newp to size of p, all // nodes must be part of the convex hull if (finalp.size() != p.size()) { if (OutLevel() > 8) std::cout << "MOERTEL: ***WRN*** MOERTEL::Overlap::ConvexHull:\n" << "MOERTEL: ***WRN*** size of convex hull " << finalp.size() << " not # nodes " << p.size() << std::endl << "MOERTEL: ***WRN*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; } // copy the polygon over to the input map p p.clear(); upper.clear(); lower.clear(); CopyPointPolygon(finalp, p); finalp.clear(); return true; }
/*----------------------------------------------------------------------* | mwgee 08/05| *----------------------------------------------------------------------*/ bool MOERTEL::Projector::ProjectNodetoSegment_Orthogonal_to_Slave( MOERTEL::Node& snode, MOERTEL::Segment& seg, double xi[], double &gap, MOERTEL::Segment& sseg) { #if 0 std::cout << "----- Projector: Node " << snode.Id() << " Segment " << seg.Id() << std::endl; std::cout << " orthogonal to Slave Segment " << sseg.Id() << std::endl; #endif if (IsTwoDimensional()) { // get the local node id of snode on sseg int lid = sseg.GetLocalNodeId(snode.Id()); if (lid<0) { std::stringstream oss; oss << "***ERR*** MOERTEL::Projector::ProjectNodetoSegment_Orthogonal_to_Slave:\n" << "***ERR*** local node id could not be found\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw ReportError(oss); } // get coordinate of local node lid double lidxi; sseg.LocalCoordinatesOfNode(lid,&lidxi); // build tangent vector in xi double g[2]; sseg.Metric(&lidxi,g,NULL); // we do a newton iteration for the projection coordinates xi // set starting value to the middle of the segment double eta = 0.0; int i = 0; double F=0.0,dF=0.0,deta=0.0; for (i=0; i<10; ++i) { F = evaluate_F_2D_SegmentOrthogonal_to_g(snode,seg,eta,gap,g); if (abs(F) < 1.0e-10) break; dF = evaluate_gradF_2D_SegmentOrthogonal_to_g(snode,seg,eta,g); deta = (-F)/dF; eta += deta; } if (abs(F)>1.0e-10) { if (OutLevel()>3) std::cout << "MOERTEL: ***WRN*** MOERTEL::Projector::ProjectNodetoSegment_Orthogonal_to_Slave:\n" << "MOERTEL: ***WRN*** Newton iteration failed to converge\n" << "MOERTEL: ***WRN*** #iterations = " << i << std::endl << "MOERTEL: ***WRN*** F(eta) = " << F << " gradF(eta) = " << dF << " eta = " << eta << " delta(eta) = " << deta << "\n" << "MOERTEL: ***WRN*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; } #if 0 std::cout << "#iterations = " << i << " F = " << F << " eta = " << eta << std::endl; #endif xi[0] = eta; return true; } else { std::stringstream oss; oss << "***ERR*** MOERTEL::Projector::ProjectNodetoSegment_Orthogonal_to_Slave:\n" << "***ERR*** 3D projection not impl.\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw ReportError(oss); } return true; }
/*----------------------------------------------------------------------* | make a triangulization of a polygon (private) mwgee 10/05| *----------------------------------------------------------------------*/ bool MOERTEL::Overlap::Triangulation() { // we have a polygon that is in clockwise order at this moment // if more than 3 points, find a center point int np = SizePointPolygon(); if (np<3) { std::cout << "***ERR*** MOERTEL::Overlap::Triangulization:\n" << "***ERR*** # point in polygon < 3 ... very strange!!!\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; exit(EXIT_FAILURE); } if (np>3) // we have to add a center point { double xi[2]; std::vector<Teuchos::RCP<MOERTEL::Point> > points; PointView(points); #if 1 // compute center point xi as centroid Centroid(xi,points,np); //std::cout << "Centroid xi : " << xi[0] << " / " << xi[1] << endl; fflush(stdout); #endif #if 0 // compute center point as nodal avergage xi[0] = xi[1] = 0.0; for (int i=0; i<np; ++i) { const double* pxi = points[i]->Xi(); xi[0] += pxi[0]; xi[1] += pxi[1]; } xi[0] /= np; xi[1] /= np; //std::cout << "Nodal xi : " << xi[0] << " / " << xi[1] << endl << endl; fflush(stdout); #endif // create a point -1 as center point and add it to the polygon AddPointtoPolygon(-1,xi); points.clear(); } // if (np>3) np = SizePointPolygon(); std::vector<Teuchos::RCP<MOERTEL::Point> > points; PointView(points); // create a MOERTEL::Node for every point int dof[3]; dof[0] = dof[1] = dof[2] = -1; // find real world coords for all points // find real world normal for all points // note that the polygon is in slave segment parameter space and is // completely contained in the slave segment. we can therefore use // slave segment values to interpolate polgon point values for (int i=0; i<np; ++i) { double x[3]; x[0] = x[1] = x[2] = 0.0; double n[3]; n[0] = n[1] = n[2] = 0.0; double val[20]; sseg_.EvaluateFunction(0,points[i]->Xi(),val,sseg_.Nnode(),NULL); MOERTEL::Node** snodes = sseg_.Nodes(); for (int j=0; j<sseg_.Nnode(); ++j) for (int k=0; k<3; ++k) { x[k] += val[j]*snodes[j]->X()[k]; n[k] += val[j]*snodes[j]->N()[k]; } const double length = MOERTEL::length(n,3); for (int j=0; j<3; ++j) n[j] /= length; // create a node with this coords and normal; MOERTEL::Node* node = new MOERTEL::Node(points[i]->Id(),x,3,dof,false,OutLevel()); node->SetN(n); // set node in point points[i]->SetNode(node); #if 0 std::cout << *points[i]; #endif } // find projection values for all points in polygon on mseg { double mxi[2]; double gap; MOERTEL::Projector projector(inter_.IsOneDimensional(),OutLevel()); for (int i=0; i<np; ++i) { Teuchos::RCP<MOERTEL::Node> node = points[i]->Node(); projector.ProjectNodetoSegment_NodalNormal(*node,mseg_,mxi,gap); // create a projected node and set it in node MOERTEL::ProjectedNode* pnode = new MOERTEL::ProjectedNode(*node,mxi,&mseg_); node->SetProjectedNode(pnode); node->SetGap(gap); #if 0 std::cout << "-------------------------------------------------------\n"; if (mseg_.Nnode()==3) { if (mxi[0]<=1. && mxi[1]<=abs(1.-mxi[0]) && mxi[0]>=0. && mxi[1]>=0.) std::cout << "OVERLAP: point " << points[i]->Id() << " is in mseg, mxi " << mxi[0] << " / " << mxi[1] << endl; else std::cout << "OVERLAP: point " << points[i]->Id() << " is NOT in mseg, mxi " << mxi[0] << " / " << mxi[1] << endl; } else if (mseg_.Nnode()==4) { if (mxi[0]<=1.001 && mxi[0]>=-1.001 && mxi[1]<=1.001 && mxi[1]>=-1.001) std::cout << "OVERLAP: point " << points[i]->Id() << " is in mseg, mxi " << mxi[0] << " / " << mxi[1] << endl; else std::cout << "OVERLAP: point " << points[i]->Id() << " is NOT in mseg, mxi " << mxi[0] << " / " << mxi[1] << endl; } else { std::cout << "***ERR*** MOERTEL::Overlap::Triangulization:\n" << "***ERR*** # nodes " << mseg_.Nnode() << " of master segment is unknown\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; exit(EXIT_FAILURE); } std::cout << "-------------------------------------------------------\n"; #endif } } // if we plan to interpolate function values at the gaussian points we need the // function values at the points if (exactvalues_==false) { // every point has 3 values of the 3 shape functions of the elements: // max 4 values from function 0 from sseg // max 4 values from function 1 from sseg // max 4 values from function 0 from mseg for (int i=0; i<np; ++i) { double val[20]; int nmval = mseg_.Nnode(); int nsval = sseg_.Nnode(); // evaluate function 0 from sseg sseg_.EvaluateFunction(0,points[i]->Xi(),val,nsval,NULL); points[i]->StoreFunctionValues(0,val,nsval); // evaluate function 1 from sseg sseg_.EvaluateFunction(1,points[i]->Xi(),val,nsval,NULL); points[i]->StoreFunctionValues(1,val,nsval); // evaluate function 0 from mseg mseg_.EvaluateFunction(0,points[i]->Node()->GetProjectedNode()->Xi(),val,nmval,NULL); points[i]->StoreFunctionValues(2,val,nmval); } } // create the triangle elements from polygon with centerpoint // In case of np==3, there is no centerpoint if (np>3) { // the polygon is in clockwise order, center point is points[0] and has // id = -1 if (points[0]->Id() != -1) { std::cout << "***ERR*** MOERTEL::Overlap::Triangulization:\n" << "***ERR*** points[0]->Id() is not -1\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; exit(EXIT_FAILURE); } int nodeid[3]; MOERTEL::Segment_BiLinearTri* tmp; MOERTEL::Function_LinearTri* func = new Function_LinearTri(OutLevel()); for (int i=2; i<np; ++i) { // there are np-1 triangles // triangle ids go from 0 to np-2 // a triangle is defined by nodes 0, i, i-1 // *this class takes ownership of triangle nodeid[0] = points[0]->Id(); nodeid[1] = points[i]->Id(); nodeid[2] = points[i-1]->Id(); tmp = new MOERTEL::Segment_BiLinearTri(i-2,3,nodeid,OutLevel()); // set a linear shape function to this triangle tmp->SetFunction(0,func); // add triangle to the *this class AddSegment(tmp->Id(),tmp); } // add the last triangle defined by nodes 0, 1, np-1 separately // *this class takes ownership of triangle nodeid[0] = points[0]->Id(); nodeid[1] = points[1]->Id(); nodeid[2] = points[np-1]->Id(); tmp = new MOERTEL::Segment_BiLinearTri(np-2,3,nodeid,OutLevel()); // set a linear shape function to this triangle tmp->SetFunction(0,func); // add triangle to the *this class AddSegment(tmp->Id(),tmp); if (func) delete func; func = NULL; } else if (np==3) // single triangle without centerpoint { int nodeid[3]; MOERTEL::Segment_BiLinearTri* tmp; MOERTEL::Function_LinearTri* func = new Function_LinearTri(OutLevel()); nodeid[0] = points[0]->Id(); nodeid[1] = points[2]->Id(); nodeid[2] = points[1]->Id(); // *this class takes ownership of triangle tmp = new MOERTEL::Segment_BiLinearTri(0,3,nodeid,OutLevel()); // set a linear shape function to this triangle tmp->SetFunction(0,func); // add triangle the *this class AddSegment(tmp->Id(),tmp); if (func) delete func; func = NULL; } else { std::cout << "***ERR*** MOERTEL::Overlap::Triangulization:\n" << "***ERR*** # point in polygon < 3 ... very strange!!!\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; exit(EXIT_FAILURE); } // create ptr topology between triangle and nodes in points std::vector<MOERTEL::Node*> nodes(np); for (int i=0; i<np; ++i) nodes[i] = points[i]->Node().get(); // loop segments and set ptr to nodes in them std::map<int,Teuchos::RCP<MOERTEL::Segment> >::iterator curr; for (curr=s_.begin(); curr != s_.end(); ++curr) curr->second->GetPtrstoNodes(nodes); nodes.clear(); points.clear(); return true; }
/*----------------------------------------------------------------------* | solve problem (public) mwgee 12/05| *----------------------------------------------------------------------*/ bool MOERTEL::Manager::Solve(Epetra_Vector& sol, const Epetra_Vector& rhs) { // test for solver parameters if (solverparams_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::Solve:\n" << "***ERR*** No solver parameters set, use SetSolverParameters(ParameterList& params)\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } // test for problemmap_ if (problemmap_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::Solve:\n" << "***ERR*** No problem map set, use SetProblemMap(Epetra_Map* map)\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } // test for inputmatrix if (inputmatrix_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::Solve:\n" << "***ERR*** No inputmatrix set, use SetInputMatrix(Epetra_CrsMatrix* inputmatrix, bool DeepCopy = false)\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } // test whether problemmap_ matches RangeMap() of inputmatrix if (!problemmap_->PointSameAs(inputmatrix_->RangeMap())) { cout << "***ERR*** MOERTEL::Manager::Solve:\n" << "***ERR*** problem map does not match range map of input matrix\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } // test whether maps of rhs and sol are ok if (!problemmap_->PointSameAs(rhs.Map()) || !problemmap_->PointSameAs(sol.Map()) ) { cout << "***ERR*** MOERTEL::Manager::Solve:\n" << "***ERR*** problem map does not match map of rhs and/or sol\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } // test whether interfaces are complete and have been integrated std::map<int,Teuchos::RCP<MOERTEL::Interface> >::iterator curr; for (curr=interface_.begin(); curr != interface_.end(); ++curr) { if (curr->second->IsComplete() == false) { cout << "***ERR*** MOERTEL::Manager::Solve:\n" << "***ERR*** interface " << curr->second->Id() << " is not IsComplete()\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } if (curr->second->IsIntegrated() == false) { cout << "***ERR*** MOERTEL::Manager::Solve:\n" << "***ERR*** interface " << curr->second->Id() << " is not integrated yet\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } } // test whether we have M and D matrix if (D_==Teuchos::null || M_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::Solve:\n" << "***ERR*** Matrix M or D is NULL\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } //--------------------------------------------------------------------------- // make solution and rhs vector matching the system Teuchos::RCP<Epetra_Vector> b = Teuchos::rcp(const_cast<Epetra_Vector*>(&rhs)); b.release(); Teuchos::RCP<Epetra_Vector> x = Teuchos::rcp(&sol); x.release(); //--------------------------------------------------------------------------- // get type of system to be used/generated Teuchos::RCP<Epetra_CrsMatrix> matrix = Teuchos::null; string system = solverparams_->get("System","None"); if (system=="None") { cout << "***WRN*** MOERTEL::Manager::Solve:\n" << "***WRN*** parameter 'System' was not chosen, using default\n" << "***WRN*** which is 'SaddleSystem'\n" << "***WRN*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; solverparams_->set("System","SaddleSystem"); system = "SaddleSystem"; } //--------------------------------------------------------------------------- // build a saddle point system if (system=="SaddleSystem" || system=="saddlesystem" || system=="SADDLESYSTEM" || system=="Saddle_System" || system=="saddle_system" || system=="SADDLE_SYSTEM") { if (saddlematrix_==Teuchos::null) { Epetra_CrsMatrix* tmp = MakeSaddleProblem(); if (!tmp) { cout << "***ERR*** MOERTEL::Manager::Solve:\n" << "***ERR*** MakeSaddleProblem() returned NULL\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } } matrix = saddlematrix_; b = Teuchos::rcp(new Epetra_Vector(*saddlemap_,true)); b.set_has_ownership(); x = Teuchos::rcp(new Epetra_Vector(*saddlemap_,false)); x.set_has_ownership(); for (int i=0; i<rhs.MyLength(); ++i) { (*b)[i] = rhs[i]; (*x)[i] = sol[i]; } } //--------------------------------------------------------------------------- // build a spd system else if (system=="SPDSystem" || system=="spdsystem" || system=="spd_system" || system=="SPD_System" || system=="SPDSYSTEM" || system=="SPD_SYSTEM") { if (spdmatrix_==Teuchos::null) { Epetra_CrsMatrix* tmp = MakeSPDProblem(); if (!tmp) { cout << "***ERR*** MOERTEL::Manager::Solve:\n" << "***ERR*** MakeSPDProblem() returned NULL\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } } matrix = spdmatrix_; // we have to multiply the rhs vector b with spdrhs_ to fit with spdmatrix_ if (spdrhs_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::Solve:\n" << "***ERR*** Cannot build righthandside for spd problem\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } Epetra_Vector *tmp = new Epetra_Vector(b->Map(),false); int err = spdrhs_->Multiply(false,*b,*tmp); if (err) { cout << "***ERR*** MOERTEL::Manager::Solve:\n" << "***ERR*** spdrhs_->Multiply returned err = " << err << endl << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } b = Teuchos::rcp(tmp); b.set_has_ownership(); } //--------------------------------------------------------------------------- // unknown parameter "System" else { cout << "***ERR*** MOERTEL::Manager::Solve:\n" << "***ERR*** Unknown type of parameter 'System': " << system << "\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return false; } //--------------------------------------------------------------------------- // create a mortar solver class instance if (solver_==Teuchos::null) solver_ = Teuchos::rcp(new MOERTEL::Solver(Comm(),OutLevel())); //--------------------------------------------------------------------------- // solve bool ok = solver_->Solve(solverparams_,matrix,x,b,*this); if (!ok) { if (Comm().MyPID()==0) cout << "***WRN*** MOERTEL::Manager::Solve:\n" << "***WRN*** MOERTEL::Solver::Solve returned an error\n" << "***WRN*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; } //--------------------------------------------------------------------------- // copy solution back to sol if neccessary if (x.has_ownership()) { for (int i=0; i<sol.MyLength(); ++i) sol[i] = (*x)[i]; } return ok; }