Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Transfer_2f.C
Go to the documentation of this file.
1 /* *******************************************************************
2  * Rocstar Simulation Suite *
3  * Copyright@2015, Illinois Rocstar LLC. All rights reserved. *
4  * *
5  * Illinois Rocstar LLC *
6  * Champaign, IL *
7  * www.illinoisrocstar.com *
8  * sales@illinoisrocstar.com *
9  * *
10  * License: See LICENSE file in top level of distribution package or *
11  * http://opensource.org/licenses/NCSA *
12  *********************************************************************/
13 /* *******************************************************************
14  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
15  * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
16  * OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
17  * NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
18  * COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
19  * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
20  * ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
21  * USE OR OTHER DEALINGS WITH THE SOFTWARE. *
22  *********************************************************************/
23 // $Id: Transfer_2f.C,v 1.9 2008/12/06 08:43:29 mtcampbe Exp $
24 
25 //===================================================================
26 // This file contains the implementation for transfering data from nodes
27 // to faces and from faces to faces.
28 //
29 // Author: Xiangmin Jiao
30 // Revision: June 15, 2001
31 //===================================================================
32 
33 #include "Transfer_2f.h"
34 #include <limits>
35 #define QUIET_NAN std::numeric_limits<Real>::quiet_NaN()
36 
38 
39 // This integrator can integrate over subset of an element.
40 // Note that the integral is added to the current values of v and area.
41 template < class _Data>
42 void
44  RFC_Pane_transfer *p_trg,
45  const _Data &data_s,
46  ENE &ene_src,
47  ENE &ene_trg,
48  int sfid_src,
49  int sfid_trg,
50  const Real alpha,
51  Facial_data &tDF,
52  Facial_data &tBF,
53  int doa)
54 {
55  // Initialize natural cooredinates of the subnodes of the subface in its
56  // parent source and target faces. Compute ncs_s only if alpha!=1.
57  Point_3 ps_s[Generic_element::MAX_SIZE], ps_t[Generic_element::MAX_SIZE];
58  Point_2 ncs_s[Generic_element::MAX_SIZE], ncs_t[Generic_element::MAX_SIZE];
59 
60  Generic_element e_s(ene_src.pane()?ene_src.size_of_edges():3,
61  ene_src.pane()?ene_src.size_of_nodes():3);
62 
63  // Initialize physical and natural coordinates for source element
64  if ( is_nodal(data_s) || alpha!=1) {
65  for ( int i=0; i<3; ++i)
66  p_src->get_nat_coor_in_element( sfid_src, i, ncs_s[i]);
67 
68  if ( alpha!=1.) {
70  Element_coor_const pnts_s( nc, p_src->coordinates(), ene_src);
71 
72  for ( int i=0; i<3; ++i) e_s.interpolate( pnts_s, ncs_s[i], &ps_s[i]);
73  }
74  }
75  else {
76  std::fill_n( &ncs_s[0][0], 6, QUIET_NAN);
77  std::fill_n( &ps_s[0][0], Generic_element::MAX_SIZE*3, QUIET_NAN);
78  }
79 
80  { // Initialize physical and natural coordinates for target element
81  for ( int i=0; i<3; ++i)
82  p_trg->get_nat_coor_in_element( sfid_trg, i, ncs_t[i]);
83 
85  Generic_element e_t( ene_trg.size_of_edges(), ene_trg.size_of_nodes());
86  Element_coor_const pnts_t( nc, p_trg->coordinates(), ene_trg);
87 
88  for ( int i=0; i<3; ++i) e_t.interpolate( pnts_t, ncs_t[i], &ps_t[i]);
89  }
90 
91  // Loop throught the quadrature points of the subfacet.
92  Vector_n vt(data_s.dimension(), 0);
93  Array_n t( vt.begin(), vt.end()), v=tDF.get_value( p_trg, ene_trg.id());
94  Real &area = tBF.get_value( p_trg, ene_trg.id())[0];
95 
96  Generic_element sub_e( 3);
97  Point_2 sub_nc, nc_s;
98 
99  if ( is_nodal( data_s)) {
100  if ( doa == 0) doa = 2; // Set the default degree of accuracy
101  }
102  else
103  doa=1;
104 
105  for ( int i=0, n = sub_e.get_num_gp(doa); i<n; ++i) {
106  sub_e.get_gp_nat_coor( i, sub_nc, doa);
107 
108  if ( is_nodal( data_s))
109  sub_e.interpolate( ncs_s, sub_nc, &nc_s);
110 
111  interpolate( e_s, data_s, nc_s, t);
112 
113  Real a=sub_e.get_gp_weight(i, doa);
114  a *= sub_e.Jacobian_det( ps_s, ps_t, alpha, sub_nc);
115 
116  t *= a;
117  v += t;
118  area += a;
119  }
120 }
121 
122 template < class _SDF>
123 void
124 Transfer_base::transfer_2f( const _SDF &sDF,
125  Facial_data &tDF,
126  const Real alpha,
127  int doa,
128  bool verbose)
129 {
130  double t0=0.;
131  if ( verbose) {
132  trg.barrier(); t0=get_wtime();
133  }
134 
135  src.replicate_data( sDF, alpha!=1);
136  // First, create buffer space for the target window and initialize
137  // the entries of the target mesh to zero.
138  trg.init_facial_buffers( tDF, 1);
139  Facial_data tBF ( trg.facial_buffer( 0));
140 
141  for ( Pane_iterator pit=trg_ps.begin(); pit!=trg_ps.end(); ++pit) {
142  Real *trg_data = (*pit)->pointer( tDF.id());
143  Real *trg_buf = (*pit)->pointer( tBF.id());
144 
145  std::fill( trg_data, trg_data+(*pit)->size_of_faces()*tDF.dimension(), 0);
146  std::fill( trg_buf, trg_buf+(*pit)->size_of_faces(), 0);
147  }
148 
149  // Second, compute the integral over the target meshes by looping through
150  // the subfaces of the target window
151  ENE ene_src, ene_trg;
152  const RFC_Pane_transfer *p_src = NULL;
153  for ( Pane_iterator pit=trg_ps.begin(); pit!=trg_ps.end(); ++pit) {
154  // Loop through the subfaces of the target window
155  for ( int i=1, size=(*pit)->size_of_subfaces(); i<=size; ++i) {
156  (*pit)->get_host_element_of_subface( i, ene_trg);
157  if ( !(*pit)->need_recv( ene_trg.id())) continue;
158 
159  const Face_ID &fid = (*pit)->get_subface_counterpart(i);
160  if ( !p_src || p_src->id()!=fid.pane_id)
161  p_src = get_src_pane( fid.pane_id);
162  if ( alpha!=1 || is_nodal( sDF.tag()) )
163  p_src->get_host_element_of_subface( fid.face_id, ene_src);
164 
165  if ( is_nodal( sDF.tag()))
166  integrate_subface( p_src, *pit,
167  make_field(sDF, p_src, ene_src),
168  ene_src, ene_trg, fid.face_id, i, alpha, tDF, tBF, doa);
169  else {
170  int id = p_src->get_parent_face( fid.face_id);
171  integrate_subface( p_src, *pit,
172  make_field(sDF, p_src, id),
173  ene_src, ene_trg, fid.face_id, i, alpha, tDF, tBF, doa);
174  }
175  }
176  }
177 
178  // Loop through the panes of the target mesh
179  for ( Pane_iterator pit=trg_ps.begin(); pit!=trg_ps.end(); ++pit) {
180  Real *trg_data = (*pit)->pointer( tDF.id());
181  Real *trg_buf = (*pit)->pointer( tBF.id());
182 
183  for ( int i=1, size=(*pit)->size_of_faces(); i<=size; ++i) {
184  if ( !(*pit)->need_recv(i)) continue;
185  tDF.get_value( trg_data, i) /= tBF.get_value( trg_buf, i)[0];
186  }
187  }
188 
189  // Clean up the transfer buffers
192 
193  if ( verbose) {
194  trg.barrier();
195  if ( trg.is_root()) {
196  std::cout << "RFACE: Transfer to faces done in "
197  << get_wtime()-t0 << " seconds." << std::endl;
198  }
199  }
200 }
201 
202 void Transfer_n2f::
204  const Real alpha, int doa, bool verb) {
205  Base::transfer_2f( sv, tf, alpha, doa, verb);
206 }
207 
208 void Transfer_f2f::
210  const Real alpha, int doa, bool verb) {
211  Base::transfer_2f( sf, tf, alpha, doa, verb);
212 }
213 
214 
216 
217 
218 
219 
220 
221 
void replicate_data(const Facial_data_const &data, bool replicate_coor)
Replicate the given data from remote processes onto local process.
void transfer(const Nodal_data_const &sv, Facial_data &tf, const Real alpha, int doa=0, bool verb=false)
The main entry to the data transfer algorithm.
Definition: Transfer_2f.C:203
An adaptor for enumerating node IDs of an element.
const RFC_Pane_transfer * get_src_pane(int i)
int face_id
the local id within the pane starting from 1.
Size dimension() const
const Real * coordinates() const
Element_var_const make_field(const Nodal_data_const &d, const RFC_Pane_transfer *pn, const ENE &ene)
Construct a element-wise accessor from nodal data and pointers.
std::vector< RFC_Pane_transfer * > trg_ps
Adpator for element-wise data container.
void integrate_subface(const RFC_Pane_transfer *p_src, RFC_Pane_transfer *p_trg, const _SDF &sDF, ENE &ene_src, ENE &ene_trg, int sfid_src, int sfid_trg, const Real alpha, Facial_data &tDF, Facial_data &tBF, int doa)
#define QUIET_NAN
Definition: Transfer_2f.C:35
int pane_id
the id of the owner pane.
double Real
Definition: mapbasic.h:322
void barrier() const
Block until all processes of the window have reached here.
SURF::Generic_element_2 Generic_element
Definition: rfc_basic.h:46
Facial_data facial_buffer(int)
int id() const
Get the local id of the element within the pane.
*********************************************************************Illinois Open Source License ****University of Illinois NCSA **Open Source License University of Illinois All rights reserved ****Developed free of to any person **obtaining a copy of this software and associated documentation to deal with the Software without including without limitation the rights to and or **sell copies of the and to permit persons to whom the **Software is furnished to do subject to the following this list of conditions and the following disclaimers ****Redistributions in binary form must reproduce the above **copyright this list of conditions and the following **disclaimers in the documentation and or other materials **provided with the distribution ****Neither the names of the Center for Simulation of Advanced the University of nor the names of its **contributors may be used to endorse or promote products derived **from this Software without specific prior written permission ****THE SOFTWARE IS PROVIDED AS WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **ARISING OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE v
Definition: roccomf90.h:20
#define RFC_END_NAME_SPACE
Definition: rfc_basic.h:29
void clear_replicated_data()
Clear all the replicate data but keep metadata.
int get_parent_face(int id) const
Get the local parent face id of a given subface.
int size_of_edges() const
Number of edges per element.
RFC_Window_transfer & trg
blockLoc i
Definition: read.cpp:79
#define RFC_BEGIN_NAME_SPACE
Definition: rfc_basic.h:28
void get_nat_coor_in_element(const int eid, const int lid, Point_2 &nc) const
Take a subface id and a local subnode id, return the natual coordinates of the subnode within the par...
void transfer_2f(const _SDF &sDF, Facial_data &tDF, const Real alpha, int doa, bool verb)
template function for transfering from nodes/faces to faces.
Definition: Transfer_2f.C:124
Value_nonconst get_value(RFC_Pane_transfer *p, int v)
int size_of_nodes() const
Number of nodes per element.
const NT & n
bool is_nodal(Tag_nodal) const
void interpolate(const Generic_element &e, const Element_var_const values, const Generic_element::Nat_coor &nc, _Value &v)
void get_host_element_of_subface(int i, Element_node_enumerator &ene) const
void init_facial_buffers(const Facial_data &nd, int)
An const adaptor for accessing nodal coordinates of a pane.
RFC_Window_transfer & src
int id() const
RFC_BEGIN_NAME_SPACE double get_wtime()
Definition: Timing.h:33
bool is_root() const
Check whether the process has rank 0.
void transfer(const Facial_data_const &sf, Facial_data &tf, const Real alpha, int doa=0, bool verb=false)
The main entry to the data transfer algorithm.
Definition: Transfer_2f.C:209
const Pane * pane() const
A global ID of a face.
std::vector< RFC_Pane_transfer * >::iterator Pane_iterator
Definition: Transfer_base.h:54